source: flexpart.git/src/concoutput_mpi.f90 @ 3481cc1

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

move license from headers to a different file

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