source: flexpart.git/src/concoutput_mpi.f90 @ 20963b1

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

Backwards deposition for the MPI version

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