source: flexpart.git/src/concoutput_surf_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@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 22.2 KB
Line 
1subroutine concoutput_surf(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  !                                                                            *
36  ! Variables:                                                                 *
37  ! outnum          number of samples                                          *
38  ! ncells          number of cells with non-zero concentrations               *
39  ! sparse          .true. if in sparse matrix format, else .false.            *
40  ! tot_mu          1 for forward, initial mass mixing ration for backw. runs  *
41  !                                                                            *
42  !*****************************************************************************
43
44  use unc_mod
45  use point_mod
46  use outg_mod
47  use par_mod
48  use com_mod
49  use mpi_mod
50  use mean_mod
51
52  implicit none
53
54  real(kind=dp) :: jul
55  integer :: itime,i,ix,jy,kz,ks,kp,l,iix,jjy,kzz,nage,jjjjmmdd,ihmmss
56  integer :: sp_count_i,sp_count_r
57  real :: sp_fact
58  real :: outnum,densityoutrecept(maxreceptor),xl,yl
59
60  !real densityoutgrid(0:numxgrid-1,0:numygrid-1,numzgrid),
61  !    +grid(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,maxpointspec_act,
62  !    +    maxageclass)
63  !real wetgrid(0:numxgrid-1,0:numygrid-1,maxspec,maxpointspec_act,
64  !    +       maxageclass)
65  !real drygrid(0:numxgrid-1,0:numygrid-1,maxspec,
66  !    +       maxpointspec_act,maxageclass)
67  !real gridsigma(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,
68  !    +       maxpointspec_act,maxageclass),
69  !    +     drygridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
70  !    +     maxpointspec_act,maxageclass),
71  !    +     wetgridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
72  !    +     maxpointspec_act,maxageclass)
73  !real factor(0:numxgrid-1,0:numygrid-1,numzgrid)
74  !real sparse_dump_r(numxgrid*numygrid*numzgrid)
75  !integer sparse_dump_i(numxgrid*numygrid*numzgrid)
76
77  !real sparse_dump_u(numxgrid*numygrid*numzgrid)
78  real(dep_prec) :: auxgrid(nclassunc)
79  real(sp) :: gridtotal,gridsigmatotal,gridtotalunc
80  real(dep_prec) :: wetgridtotal,wetgridsigmatotal,wetgridtotalunc
81  real(dep_prec) :: drygridtotal,drygridsigmatotal,drygridtotalunc
82  real :: halfheight,dz,dz1,dz2,tot_mu(maxspec,maxpointspec_act)
83  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
84  real,parameter :: weightair=28.97
85  logical :: sp_zer
86  character :: adate*8,atime*6
87  character(len=3) :: anspec
88  integer :: mind
89! mind        eso:added to get consistent results between 2&3-fields versions
90
91! Measure execution time
92  if (mp_measure_time) call mpif_mtime('rootonly',0)
93
94
95  ! Determine current calendar date, needed for the file name
96  !**********************************************************
97
98  jul=bdate+real(itime,kind=dp)/86400._dp
99  call caldate(jul,jjjjmmdd,ihmmss)
100  write(adate,'(i8.8)') jjjjmmdd
101  write(atime,'(i6.6)') ihmmss
102  !write(unitdates,'(a)') adate//atime
103
104    open(unitdates,file=path(2)(1:length(2))//'dates', ACCESS='APPEND')
105      write(unitdates,'(a)') adate//atime
106    close(unitdates)
107
108  ! For forward simulations, output fields have dimension MAXSPEC,
109  ! for backward simulations, output fields have dimension MAXPOINT.
110  ! Thus, make loops either about nspec, or about numpoint
111  !*****************************************************************
112
113
114  if (ldirect.eq.1) then
115    do ks=1,nspec
116      do kp=1,maxpointspec_act
117        tot_mu(ks,kp)=1
118      end do
119    end do
120  else
121    do ks=1,nspec
122      do kp=1,maxpointspec_act
123        tot_mu(ks,kp)=xmass(kp,ks)
124      end do
125    end do
126  endif
127
128
129  if (verbosity.eq.1) then
130     print*,'concoutput_surf 2'
131     CALL SYSTEM_CLOCK(count_clock)
132     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
133  endif
134
135  !*******************************************************************
136  ! Compute air density: sufficiently accurate to take it
137  ! from coarse grid at some time
138  ! Determine center altitude of output layer, and interpolate density
139  ! data to that altitude
140  !*******************************************************************
141
142  mind=memind(2)
143  do kz=1,numzgrid
144    if (kz.eq.1) then
145      halfheight=outheight(1)/2.
146    else
147      halfheight=(outheight(kz)+outheight(kz-1))/2.
148    endif
149    do kzz=2,nz
150      if ((height(kzz-1).lt.halfheight).and. &
151           (height(kzz).gt.halfheight)) goto 46
152    end do
15346   kzz=max(min(kzz,nz),2)
154    dz1=halfheight-height(kzz-1)
155    dz2=height(kzz)-halfheight
156    dz=dz1+dz2
157    do jy=0,numygrid-1
158      do ix=0,numxgrid-1
159        xl=outlon0+real(ix)*dxout
160        yl=outlat0+real(jy)*dyout
161        xl=(xl-xlon0)/dx
162        yl=(yl-ylat0)/dy
163        iix=max(min(nint(xl),nxmin1),0)
164        jjy=max(min(nint(yl),nymin1),0)
165        ! densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
166        !      rho(iix,jjy,kzz-1,2)*dz2)/dz
167        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,mind)*dz1+ &
168             rho(iix,jjy,kzz-1,mind)*dz2)/dz
169      end do
170    end do
171  end do
172
173  do i=1,numreceptor
174    xl=xreceptor(i)
175    yl=yreceptor(i)
176    iix=max(min(nint(xl),nxmin1),0)
177    jjy=max(min(nint(yl),nymin1),0)
178    !densityoutrecept(i)=rho(iix,jjy,1,2)
179    densityoutrecept(i)=rho(iix,jjy,1,mind)
180  end do
181
182
183  ! Output is different for forward and backward simulations
184    do kz=1,numzgrid
185      do jy=0,numygrid-1
186        do ix=0,numxgrid-1
187          if (ldirect.eq.1) then
188            factor3d(ix,jy,kz)=1.e12/volume(ix,jy,kz)/outnum
189          else
190            factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
191          endif
192        end do
193      end do
194    end do
195
196  !*********************************************************************
197  ! Determine the standard deviation of the mean concentration or mixing
198  ! ratio (uncertainty of the output) and the dry and wet deposition
199  !*********************************************************************
200
201  if (verbosity.eq.1) then
202     print*,'concoutput_surf 3 (sd)'
203     CALL SYSTEM_CLOCK(count_clock)
204     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
205  endif
206  gridtotal=0.
207  gridsigmatotal=0.
208  gridtotalunc=0.
209  wetgridtotal=0.
210  wetgridsigmatotal=0.
211  wetgridtotalunc=0.
212  drygridtotal=0.
213  drygridsigmatotal=0.
214  drygridtotalunc=0.
215
216  do ks=1,nspec
217
218  write(anspec,'(i3.3)') ks
219  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
220    if (ldirect.eq.1) then
221      open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_'//adate// &
222           atime//'_'//anspec,form='unformatted')
223    else
224      open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_'//adate// &
225           atime//'_'//anspec,form='unformatted')
226    endif
227     write(unitoutgrid) itime
228   endif
229
230  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
231   open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_'//adate// &
232        atime//'_'//anspec,form='unformatted')
233
234    write(unitoutgridppt) itime
235  endif
236
237  do kp=1,maxpointspec_act
238  do nage=1,nageclass
239
240    do jy=0,numygrid-1
241      do ix=0,numxgrid-1
242
243  ! WET DEPOSITION
244        if ((WETDEP).and.(ldirect.gt.0)) then
245            do l=1,nclassunc
246              auxgrid(l)=wetgridunc0(ix,jy,ks,kp,l,nage)
247            end do
248            call mean(auxgrid,wetgrid(ix,jy), &
249                 wetgridsigma(ix,jy),nclassunc)
250  ! Multiply by number of classes to get total concentration
251            wetgrid(ix,jy)=wetgrid(ix,jy) &
252                 *nclassunc
253            wetgridtotal=wetgridtotal+wetgrid(ix,jy)
254  ! Calculate standard deviation of the mean
255            wetgridsigma(ix,jy)= &
256                 wetgridsigma(ix,jy)* &
257                 sqrt(real(nclassunc))
258            wetgridsigmatotal=wetgridsigmatotal+ &
259                 wetgridsigma(ix,jy)
260        endif
261
262  ! DRY DEPOSITION
263        if ((DRYDEP).and.(ldirect.gt.0)) then
264            do l=1,nclassunc
265              auxgrid(l)=drygridunc0(ix,jy,ks,kp,l,nage)
266            end do
267            call mean(auxgrid,drygrid(ix,jy), &
268                 drygridsigma(ix,jy),nclassunc)
269  ! Multiply by number of classes to get total concentration
270            drygrid(ix,jy)=drygrid(ix,jy)* &
271                 nclassunc
272            drygridtotal=drygridtotal+drygrid(ix,jy)
273  ! Calculate standard deviation of the mean
274            drygridsigma(ix,jy)= &
275                 drygridsigma(ix,jy)* &
276                 sqrt(real(nclassunc))
277125         drygridsigmatotal=drygridsigmatotal+ &
278                 drygridsigma(ix,jy)
279        endif
280
281  ! CONCENTRATION OR MIXING RATIO
282        do kz=1,numzgrid
283            do l=1,nclassunc
284              auxgrid(l)=gridunc(ix,jy,kz,ks,kp,l,nage)
285            end do
286            call mean(auxgrid,grid(ix,jy,kz), &
287                 gridsigma(ix,jy,kz),nclassunc)
288  ! Multiply by number of classes to get total concentration
289            grid(ix,jy,kz)= &
290                 grid(ix,jy,kz)*nclassunc
291            gridtotal=gridtotal+grid(ix,jy,kz)
292  ! Calculate standard deviation of the mean
293            gridsigma(ix,jy,kz)= &
294                 gridsigma(ix,jy,kz)* &
295                 sqrt(real(nclassunc))
296            gridsigmatotal=gridsigmatotal+ &
297                 gridsigma(ix,jy,kz)
298        end do
299      end do
300    end do
301
302
303  !*******************************************************************
304  ! Generate output: may be in concentration (ng/m3) or in mixing
305  ! ratio (ppt) or both
306  ! Output the position and the values alternated multiplied by
307  ! 1 or -1, first line is number of values, number of positions
308  ! For backward simulations, the unit is seconds, stored in grid_time
309  !*******************************************************************
310
311  if (verbosity.eq.1) then
312     print*,'concoutput_surf 4 (output)'
313     CALL SYSTEM_CLOCK(count_clock)
314     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
315  endif
316
317  ! Concentration output
318  !*********************
319
320  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
321
322  if (verbosity.eq.1) then
323     print*,'concoutput_surf (Wet deposition)'
324     CALL SYSTEM_CLOCK(count_clock)
325     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
326  endif
327
328  ! Wet deposition
329         sp_count_i=0
330         sp_count_r=0
331         sp_fact=-1.
332         sp_zer=.true.
333         if ((ldirect.eq.1).and.(WETDEP)) then
334         do jy=0,numygrid-1
335            do ix=0,numxgrid-1
336  ! concentraion greater zero
337              if (wetgrid(ix,jy).gt.smallnum) then
338                 if (sp_zer.eqv..true.) then ! first non zero value
339                    sp_count_i=sp_count_i+1
340                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
341                    sp_zer=.false.
342                    sp_fact=sp_fact*(-1.)
343                 endif
344                 sp_count_r=sp_count_r+1
345                 sparse_dump_r(sp_count_r)= &
346                      sp_fact*1.e12*wetgrid(ix,jy)/area(ix,jy)
347                 sparse_dump_u(sp_count_r)= &
348                      1.e12*wetgridsigma(ix,jy)/area(ix,jy)
349              else ! concentration is zero
350                  sp_zer=.true.
351              endif
352            end do
353         end do
354         else
355            sp_count_i=0
356            sp_count_r=0
357         endif
358         write(unitoutgrid) sp_count_i
359         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
360         write(unitoutgrid) sp_count_r
361         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
362!         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
363
364  if (verbosity.eq.1) then
365     print*,'concoutput_surf (Dry deposition)'
366     CALL SYSTEM_CLOCK(count_clock)
367     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
368  endif
369  ! Dry deposition
370         sp_count_i=0
371         sp_count_r=0
372         sp_fact=-1.
373         sp_zer=.true.
374         if ((ldirect.eq.1).and.(DRYDEP)) then
375          do jy=0,numygrid-1
376            do ix=0,numxgrid-1
377              if (drygrid(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* &
387                      1.e12*drygrid(ix,jy)/area(ix,jy)
388                  sparse_dump_u(sp_count_r)= &
389                      1.e12*drygridsigma(ix,jy)/area(ix,jy)
390              else ! concentration is zero
391                  sp_zer=.true.
392              endif
393            end do
394          end do
395         else
396            sp_count_i=0
397            sp_count_r=0
398         endif
399         write(unitoutgrid) sp_count_i
400         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
401         write(unitoutgrid) sp_count_r
402         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
403!         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
404
405  if (verbosity.eq.1) then
406     print*,'concoutput_surf (Concentrations)'
407     CALL SYSTEM_CLOCK(count_clock)
408     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
409  endif
410
411  ! Concentrations
412
413  ! surf_only write only 1st layer
414
415         sp_count_i=0
416         sp_count_r=0
417         sp_fact=-1.
418         sp_zer=.true.
419          do kz=1,1
420            do jy=0,numygrid-1
421              do ix=0,numxgrid-1
422                if (grid(ix,jy,kz).gt.smallnum) then
423                  if (sp_zer.eqv..true.) then ! first non zero value
424                    sp_count_i=sp_count_i+1
425                    sparse_dump_i(sp_count_i)= &
426                         ix+jy*numxgrid+kz*numxgrid*numygrid
427                    sp_zer=.false.
428                    sp_fact=sp_fact*(-1.)
429                  endif
430                   sp_count_r=sp_count_r+1
431                   sparse_dump_r(sp_count_r)= &
432                        sp_fact* &
433                        grid(ix,jy,kz)* &
434                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
435  !                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
436  !    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
437                   sparse_dump_u(sp_count_r)= &
438                        gridsigma(ix,jy,kz)* &
439                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
440                else ! concentration is zero
441                  sp_zer=.true.
442                endif
443             end do
444           end do
445         end do
446         write(unitoutgrid) sp_count_i
447         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
448         write(unitoutgrid) sp_count_r
449         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
450!         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
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)/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) (sparse_dump_u(i),i=1,sp_count_r)
495
496
497  ! Dry deposition
498         sp_count_i=0
499         sp_count_r=0
500         sp_fact=-1.
501         sp_zer=.true.
502         if ((ldirect.eq.1).and.(DRYDEP)) then
503          do jy=0,numygrid-1
504            do ix=0,numxgrid-1
505                if (drygrid(ix,jy).gt.smallnum) then
506                  if (sp_zer.eqv..true.) then ! first non zero value
507                    sp_count_i=sp_count_i+1
508                    sparse_dump_i(sp_count_i)= &
509                         ix+jy*numxgrid
510                    sp_zer=.false.
511                    sp_fact=sp_fact*(-1)
512                 endif
513                 sp_count_r=sp_count_r+1
514                 sparse_dump_r(sp_count_r)= &
515                      sp_fact* &
516                      1.e12*drygrid(ix,jy)/area(ix,jy)
517                 sparse_dump_u(sp_count_r)= &
518                      1.e12*drygridsigma(ix,jy)/area(ix,jy)
519              else ! concentration is zero
520                  sp_zer=.true.
521              endif
522            end do
523          end do
524         else
525           sp_count_i=0
526           sp_count_r=0
527         endif
528         write(unitoutgridppt) sp_count_i
529         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
530         write(unitoutgridppt) sp_count_r
531         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
532!         write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
533
534
535  ! Mixing ratios
536
537  ! surf_only write only 1st layer
538
539         sp_count_i=0
540         sp_count_r=0
541         sp_fact=-1.
542         sp_zer=.true.
543          do kz=1,1
544            do jy=0,numygrid-1
545              do ix=0,numxgrid-1
546                if (grid(ix,jy,kz).gt.smallnum) then
547                  if (sp_zer.eqv..true.) then ! first non zero value
548                    sp_count_i=sp_count_i+1
549                    sparse_dump_i(sp_count_i)= &
550                         ix+jy*numxgrid+kz*numxgrid*numygrid
551                    sp_zer=.false.
552                    sp_fact=sp_fact*(-1.)
553                 endif
554                 sp_count_r=sp_count_r+1
555                 sparse_dump_r(sp_count_r)= &
556                      sp_fact* &
557                      1.e12*grid(ix,jy,kz) &
558                      /volume(ix,jy,kz)/outnum* &
559                      weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
560                 sparse_dump_u(sp_count_r)= &
561                      1.e12*gridsigma(ix,jy,kz)/volume(ix,jy,kz)/ &
562                      outnum*weightair/weightmolar(ks)/ &
563                      densityoutgrid(ix,jy,kz)
564              else ! concentration is zero
565                  sp_zer=.true.
566              endif
567              end do
568            end do
569          end do
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) (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_surf
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG