source: flexpart.git/src/concoutput_mpi.f90 @ 478e9e6

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

dates file is now created at program start, instead of appending to it if it already exists

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