source: flexpart.git/src/concoutput_mpi.f90

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file was 92fab65, checked in by Ignacio Pisso <ip@…>, 4 years ago

add SPDX-License-Identifier to all .f90 files

  • Property mode set to 100644
File size: 23.0 KB
Line 
1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
3
4subroutine concoutput(itime,outnum,gridtotalunc,wetgridtotalunc, &
5     drygridtotalunc)
6!                        i     i          o             o
7!       o
8!*****************************************************************************
9!                                                                            *
10!     Output of the concentration grid and the receptor concentrations.      *
11!                                                                            *
12!     Author: A. Stohl                                                       *
13!                                                                            *
14!     24 May 1995                                                            *
15!                                                                            *
16!     13 April 1999, Major update: if output size is smaller, dump output    *
17!                    in sparse matrix format; additional output of           *
18!                    uncertainty                                             *
19!                                                                            *
20!     05 April 2000, Major update: output of age classes; output for backward*
21!                    runs is time spent in grid cell times total mass of     *
22!                    species.                                                *
23!                                                                            *
24!     17 February 2002, Appropriate dimensions for backward and forward runs *
25!                       are now specified in file par_mod                    *
26!                                                                            *
27!     June 2006, write grid in sparse matrix with a single write command     *
28!                in order to save disk space                                 *
29!                                                                            *
30!     2008 new sparse matrix format                                          *
31!                                                                            *
32!     Changes eso:                                                           *
33!     2014        MPI version    This routine is only called by root MPI     *
34!                                process (the other processes have sent      *
35!                                their fields to root)                       *
36!*****************************************************************************
37!                                                                            *
38! Variables:                                                                 *
39! outnum          number of samples                                          *
40! ncells          number of cells with non-zero concentrations               *
41! sparse          .true. if in sparse matrix format, else .false.            *
42! tot_mu          1 for forward, initial mass mixing ration for backw. runs  *
43!                                                                            *
44!*****************************************************************************
45
46  use unc_mod
47  use point_mod
48  use outg_mod
49  use par_mod
50  use com_mod
51  use mpi_mod
52  use mean_mod
53
54  implicit none
55
56  real(kind=dp) :: jul
57  integer :: itime,i,ix,jy,kz,ks,kp,l,iix,jjy,kzz,nage,jjjjmmdd,ihmmss
58  integer :: sp_count_i,sp_count_r
59  real :: sp_fact
60  real :: outnum,densityoutrecept(maxreceptor),xl,yl
61
62!real densityoutgrid(0:numxgrid-1,0:numygrid-1,numzgrid),
63!    +grid(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,maxpointspec_act,
64!    +    maxageclass)
65!real wetgrid(0:numxgrid-1,0:numygrid-1,maxspec,maxpointspec_act,
66!    +       maxageclass)
67!real drygrid(0:numxgrid-1,0:numygrid-1,maxspec,
68!    +       maxpointspec_act,maxageclass)
69!real gridsigma(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,
70!    +       maxpointspec_act,maxageclass),
71!    +     drygridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
72!    +     maxpointspec_act,maxageclass),
73!    +     wetgridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
74!    +     maxpointspec_act,maxageclass)
75!real factor(0:numxgrid-1,0:numygrid-1,numzgrid)
76!real sparse_dump_r(numxgrid*numygrid*numzgrid)
77!integer sparse_dump_i(numxgrid*numygrid*numzgrid)
78
79!real sparse_dump_u(numxgrid*numygrid*numzgrid)
80  real(dep_prec) :: auxgrid(nclassunc)
81  real(sp) :: gridtotal,gridsigmatotal,gridtotalunc
82  real(dep_prec) :: wetgridtotal,wetgridsigmatotal,wetgridtotalunc
83  real(dep_prec) :: drygridtotal,drygridsigmatotal,drygridtotalunc
84  real :: halfheight,dz,dz1,dz2,tot_mu(maxspec,maxpointspec_act)
85  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
86  real,parameter :: weightair=28.97
87  logical :: sp_zer
88  logical,save :: init=.true.
89  character :: adate*8,atime*6
90  character(len=3) :: anspec
91  integer :: mind
92! mind        eso: added to ensure identical results between 2&3-fields versions
93  character(LEN=8),save :: file_stat='REPLACE'
94  logical :: ldates_file
95  integer :: ierr
96  character(LEN=100) :: dates_char
97
98! Measure execution time
99  if (mp_measure_time) call mpif_mtime('rootonly',0)
100
101! Determine current calendar date, needed for the file name
102!**********************************************************
103
104  jul=bdate+real(itime,kind=dp)/86400._dp
105  call caldate(jul,jjjjmmdd,ihmmss)
106  write(adate,'(i8.8)') jjjjmmdd
107  write(atime,'(i6.6)') ihmmss
108
109! Overwrite existing dates file on first call, later append to it
110! This fixes a bug where the dates file kept growing across multiple runs
111
112! If 'dates' file exists in output directory, make a backup
113  inquire(file=path(2)(1:length(2))//'dates', exist=ldates_file)
114  if (ldates_file.and.init) then
115    open(unit=unitdates, file=path(2)(1:length(2))//'dates',form='formatted', &
116         &access='sequential', status='old', action='read', iostat=ierr)
117    open(unit=unittmp, file=path(2)(1:length(2))//'dates.bak', access='sequential', &
118         &status='replace', action='write', form='formatted', iostat=ierr)
119    do while (.true.)
120      read(unitdates, '(a)', iostat=ierr) dates_char
121      if (ierr.ne.0) exit
122      write(unit=unittmp, fmt='(a)', iostat=ierr, advance='yes') trim(dates_char)
123!      if (ierr.ne.0) write(*,*) "Write error, ", ierr
124    end do
125    close(unit=unitdates)
126    close(unit=unittmp)
127  end if
128
129  open(unitdates,file=path(2)(1:length(2))//'dates', ACCESS='APPEND', STATUS=file_stat)
130  write(unitdates,'(a)') adate//atime
131  close(unitdates) 
132
133  ! Overwrite existing dates file on first call, later append to it
134  ! This fixes a bug where the dates file kept growing across multiple runs
135  IF (init) THEN
136    file_stat='OLD'
137    init=.false.
138  END IF
139
140
141! For forward simulations, output fields have dimension MAXSPEC,
142! for backward simulations, output fields have dimension MAXPOINT.
143! Thus, make loops either about nspec, or about numpoint
144!*****************************************************************
145
146
147  if (ldirect.eq.1) then
148    do ks=1,nspec
149      do kp=1,maxpointspec_act
150        tot_mu(ks,kp)=1
151      end do
152    end do
153  else
154    do ks=1,nspec
155      do kp=1,maxpointspec_act
156        tot_mu(ks,kp)=xmass(kp,ks)
157      end do
158    end do
159  endif
160
161
162!*******************************************************************
163! Compute air density: sufficiently accurate to take it
164! from coarse grid at some time
165! Determine center altitude of output layer, and interpolate density
166! data to that altitude
167!*******************************************************************
168
169  mind=memind(2)
170  do kz=1,numzgrid
171    if (kz.eq.1) then
172      halfheight=outheight(1)/2.
173    else
174      halfheight=(outheight(kz)+outheight(kz-1))/2.
175    endif
176    do kzz=2,nz
177      if ((height(kzz-1).lt.halfheight).and. &
178           (height(kzz).gt.halfheight)) goto 46
179    end do
18046  kzz=max(min(kzz,nz),2)
181    dz1=halfheight-height(kzz-1)
182    dz2=height(kzz)-halfheight
183    dz=dz1+dz2
184    do jy=0,numygrid-1
185      do ix=0,numxgrid-1
186        xl=outlon0+real(ix)*dxout
187        yl=outlat0+real(jy)*dyout
188        xl=(xl-xlon0)/dx
189        yl=(yl-ylat0)/dy !v9.1.1
190        iix=max(min(nint(xl),nxmin1),0)
191        jjy=max(min(nint(yl),nymin1),0)
192        ! densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
193        !      rho(iix,jjy,kzz-1,2)*dz2)/dz
194        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,mind)*dz1+ &
195             rho(iix,jjy,kzz-1,mind)*dz2)/dz
196      end do
197    end do
198  end do
199
200  do i=1,numreceptor
201    xl=xreceptor(i)
202    yl=yreceptor(i)
203    iix=max(min(nint(xl),nxmin1),0)
204    jjy=max(min(nint(yl),nymin1),0)
205    !densityoutrecept(i)=rho(iix,jjy,1,2)
206    densityoutrecept(i)=rho(iix,jjy,1,mind)
207  end do
208
209
210! Output is different for forward and backward simulations
211  do kz=1,numzgrid
212    do jy=0,numygrid-1
213      do ix=0,numxgrid-1
214        if (ldirect.eq.1) then
215          factor3d(ix,jy,kz)=1.e12/volume(ix,jy,kz)/outnum
216        else
217          factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
218        endif
219      end do
220    end do
221  end do
222
223!*********************************************************************
224! Determine the standard deviation of the mean concentration or mixing
225! ratio (uncertainty of the output) and the dry and wet deposition
226!*********************************************************************
227
228  gridtotal=0.
229  gridsigmatotal=0.
230  gridtotalunc=0.
231  wetgridtotal=0.
232  wetgridsigmatotal=0.
233  wetgridtotalunc=0.
234  drygridtotal=0.
235  drygridsigmatotal=0.
236  drygridtotalunc=0.
237
238  do ks=1,nspec
239
240    write(anspec,'(i3.3)') ks
241
242    if (DRYBKDEP.or.WETBKDEP) then !scavdep output
243      if (DRYBKDEP) &
244      open(unitoutgrid,file=path(2)(1:length(2))//'grid_drydep_'//adate// &
245           atime//'_'//anspec,form='unformatted')
246      if (WETBKDEP) &
247      open(unitoutgrid,file=path(2)(1:length(2))//'grid_wetdep_'//adate// &
248           atime//'_'//anspec,form='unformatted')
249      write(unitoutgrid) itime
250    else
251      if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
252        if (ldirect.eq.1) then
253          open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_'//adate// &
254               atime//'_'//anspec,form='unformatted')
255        else
256          open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_'//adate// &
257               atime//'_'//anspec,form='unformatted')
258        endif
259        write(unitoutgrid) itime
260      endif
261
262      if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
263        open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_'//adate// &
264             atime//'_'//anspec,form='unformatted')
265
266        write(unitoutgridppt) itime
267      endif
268    endif ! if deposition output
269
270    do kp=1,maxpointspec_act
271      do nage=1,nageclass
272
273        do jy=0,numygrid-1
274          do ix=0,numxgrid-1
275
276! WET DEPOSITION
277            if ((WETDEP).and.(ldirect.gt.0)) then
278              do l=1,nclassunc
279                auxgrid(l)=wetgridunc0(ix,jy,ks,kp,l,nage)
280              end do
281              call mean(auxgrid,wetgrid(ix,jy), &
282                   wetgridsigma(ix,jy),nclassunc)
283! Multiply by number of classes to get total concentration
284              wetgrid(ix,jy)=wetgrid(ix,jy) &
285                   *nclassunc
286              wetgridtotal=wetgridtotal+wetgrid(ix,jy)
287! Calculate standard deviation of the mean
288              wetgridsigma(ix,jy)= &
289                   wetgridsigma(ix,jy)* &
290                   sqrt(real(nclassunc))
291              wetgridsigmatotal=wetgridsigmatotal+ &
292                   wetgridsigma(ix,jy)
293            endif
294
295! DRY DEPOSITION
296            if ((DRYDEP).and.(ldirect.gt.0)) then
297              do l=1,nclassunc
298                auxgrid(l)=drygridunc0(ix,jy,ks,kp,l,nage)
299              end do
300              call mean(auxgrid,drygrid(ix,jy), &
301                   drygridsigma(ix,jy),nclassunc)
302! Multiply by number of classes to get total concentration
303              drygrid(ix,jy)=drygrid(ix,jy)* &
304                   nclassunc
305              drygridtotal=drygridtotal+drygrid(ix,jy)
306! Calculate standard deviation of the mean
307              drygridsigma(ix,jy)= &
308                   drygridsigma(ix,jy)* &
309                   sqrt(real(nclassunc))
310              drygridsigmatotal=drygridsigmatotal+ &
311                   drygridsigma(ix,jy)
312            endif
313
314! CONCENTRATION OR MIXING RATIO
315            do kz=1,numzgrid
316              do l=1,nclassunc
317                auxgrid(l)=gridunc(ix,jy,kz,ks,kp,l,nage)
318              end do
319              call mean(auxgrid,grid(ix,jy,kz), &
320                   gridsigma(ix,jy,kz),nclassunc)
321! Multiply by number of classes to get total concentration
322              grid(ix,jy,kz)= &
323                   grid(ix,jy,kz)*nclassunc
324              gridtotal=gridtotal+grid(ix,jy,kz)
325! Calculate standard deviation of the mean
326              gridsigma(ix,jy,kz)= &
327                   gridsigma(ix,jy,kz)* &
328                   sqrt(real(nclassunc))
329              gridsigmatotal=gridsigmatotal+ &
330                   gridsigma(ix,jy,kz)
331            end do
332          end do
333        end do
334
335
336
337
338!*******************************************************************
339! Generate output: may be in concentration (ng/m3) or in mixing
340! ratio (ppt) or both
341! Output the position and the values alternated multiplied by
342! 1 or -1, first line is number of values, number of positions
343! For backward simulations, the unit is seconds, stored in grid_time
344!*******************************************************************
345
346! Concentration output
347!*********************
348        if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
349
350! Wet deposition
351          sp_count_i=0
352          sp_count_r=0
353          sp_fact=-1.
354          sp_zer=.true.
355          if ((ldirect.eq.1).and.(WETDEP)) then
356            do jy=0,numygrid-1
357              do ix=0,numxgrid-1
358!oncentraion greater zero
359                if (wetgrid(ix,jy).gt.smallnum) then
360                  if (sp_zer.eqv..true.) then ! first non zero value
361                    sp_count_i=sp_count_i+1
362                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
363                    sp_zer=.false.
364                    sp_fact=sp_fact*(-1.)
365                  endif
366                  sp_count_r=sp_count_r+1
367                  sparse_dump_r(sp_count_r)= &
368                       sp_fact*1.e12*wetgrid(ix,jy)/area(ix,jy)
369!                sparse_dump_u(sp_count_r)=
370!+                1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
371                else ! concentration is zero
372                  sp_zer=.true.
373                endif
374              end do
375            end do
376          else
377            sp_count_i=0
378            sp_count_r=0
379          endif
380          write(unitoutgrid) sp_count_i
381          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
382          write(unitoutgrid) sp_count_r
383          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
384!       write(unitoutgrid) sp_count_u
385!       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
386
387! Dry deposition
388          sp_count_i=0
389          sp_count_r=0
390          sp_fact=-1.
391          sp_zer=.true.
392          if ((ldirect.eq.1).and.(DRYDEP)) then
393            do jy=0,numygrid-1
394              do ix=0,numxgrid-1
395                if (drygrid(ix,jy).gt.smallnum) then
396                  if (sp_zer.eqv..true.) then ! first non zero value
397                    sp_count_i=sp_count_i+1
398                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
399                    sp_zer=.false.
400                    sp_fact=sp_fact*(-1.)
401                  endif
402                  sp_count_r=sp_count_r+1
403                  sparse_dump_r(sp_count_r)= &
404                       sp_fact* &
405                       1.e12*drygrid(ix,jy)/area(ix,jy)
406!                sparse_dump_u(sp_count_r)=
407!+                1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
408                else ! concentration is zero
409                  sp_zer=.true.
410                endif
411              end do
412            end do
413          else
414            sp_count_i=0
415            sp_count_r=0
416          endif
417          write(unitoutgrid) sp_count_i
418          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
419          write(unitoutgrid) sp_count_r
420          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
421!       write(*,*) sp_count_u
422!       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
423
424
425
426! Concentrations
427          sp_count_i=0
428          sp_count_r=0
429          sp_fact=-1.
430          sp_zer=.true.
431          do kz=1,numzgrid
432            do jy=0,numygrid-1
433              do ix=0,numxgrid-1
434                if (grid(ix,jy,kz).gt.smallnum) then
435                  if (sp_zer.eqv..true.) then ! first non zero value
436                    sp_count_i=sp_count_i+1
437                    sparse_dump_i(sp_count_i)= &
438                         ix+jy*numxgrid+kz*numxgrid*numygrid
439                    sp_zer=.false.
440                    sp_fact=sp_fact*(-1.)
441                  endif
442                  sp_count_r=sp_count_r+1
443                  if (lparticlecountoutput) then
444                    sparse_dump_r(sp_count_r)= &
445                         sp_fact* &
446                         grid(ix,jy,kz)
447                  else
448                  sparse_dump_r(sp_count_r)= &
449                       sp_fact* &
450                       grid(ix,jy,kz)* &
451                       factor3d(ix,jy,kz)/tot_mu(ks,kp)
452                  end if
453
454
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  creceptor(:,:)=0.
641  gridunc(:,:,:,:,:,:,:)=0.
642
643  if (mp_measure_time) call mpif_mtime('rootonly',1)
644 
645end subroutine concoutput
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG