source: flexpart.git/src/concoutput_surf_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: 22.3 KB
Line 
1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
3
4subroutine concoutput_surf(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  !                                                                            *
39  ! Variables:                                                                 *
40  ! outnum          number of samples                                          *
41  ! ncells          number of cells with non-zero concentrations               *
42  ! sparse          .true. if in sparse matrix format, else .false.            *
43  ! tot_mu          1 for forward, initial mass mixing ration for backw. runs  *
44  !                                                                            *
45  !*****************************************************************************
46
47  use unc_mod
48  use point_mod
49  use outg_mod
50  use par_mod
51  use com_mod
52  use mpi_mod
53  use mean_mod
54
55  implicit none
56
57  real(kind=dp) :: jul
58  integer :: itime,i,ix,jy,kz,ks,kp,l,iix,jjy,kzz,nage,jjjjmmdd,ihmmss
59  integer :: sp_count_i,sp_count_r
60  real :: sp_fact
61  real :: outnum,densityoutrecept(maxreceptor),xl,yl
62
63  !real densityoutgrid(0:numxgrid-1,0:numygrid-1,numzgrid),
64  !    +grid(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,maxpointspec_act,
65  !    +    maxageclass)
66  !real wetgrid(0:numxgrid-1,0:numygrid-1,maxspec,maxpointspec_act,
67  !    +       maxageclass)
68  !real drygrid(0:numxgrid-1,0:numygrid-1,maxspec,
69  !    +       maxpointspec_act,maxageclass)
70  !real gridsigma(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,
71  !    +       maxpointspec_act,maxageclass),
72  !    +     drygridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
73  !    +     maxpointspec_act,maxageclass),
74  !    +     wetgridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
75  !    +     maxpointspec_act,maxageclass)
76  !real factor(0:numxgrid-1,0:numygrid-1,numzgrid)
77  !real sparse_dump_r(numxgrid*numygrid*numzgrid)
78  !integer sparse_dump_i(numxgrid*numygrid*numzgrid)
79
80  !real sparse_dump_u(numxgrid*numygrid*numzgrid)
81  real(dep_prec) :: auxgrid(nclassunc)
82  real(sp) :: gridtotal,gridsigmatotal,gridtotalunc
83  real(dep_prec) :: wetgridtotal,wetgridsigmatotal,wetgridtotalunc
84  real(dep_prec) :: drygridtotal,drygridsigmatotal,drygridtotalunc
85  real :: halfheight,dz,dz1,dz2,tot_mu(maxspec,maxpointspec_act)
86  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
87  real,parameter :: weightair=28.97
88  logical :: sp_zer
89  character :: adate*8,atime*6
90  character(len=3) :: anspec
91  integer :: mind
92! mind        eso:added to get consistent results between 2&3-fields versions
93
94! Measure execution time
95  if (mp_measure_time) call mpif_mtime('rootonly',0)
96
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
105  !write(unitdates,'(a)') adate//atime
106
107    open(unitdates,file=path(2)(1:length(2))//'dates', ACCESS='APPEND')
108      write(unitdates,'(a)') adate//atime
109    close(unitdates)
110
111  ! For forward simulations, output fields have dimension MAXSPEC,
112  ! for backward simulations, output fields have dimension MAXPOINT.
113  ! Thus, make loops either about nspec, or about numpoint
114  !*****************************************************************
115
116
117  if (ldirect.eq.1) then
118    do ks=1,nspec
119      do kp=1,maxpointspec_act
120        tot_mu(ks,kp)=1
121      end do
122    end do
123  else
124    do ks=1,nspec
125      do kp=1,maxpointspec_act
126        tot_mu(ks,kp)=xmass(kp,ks)
127      end do
128    end do
129  endif
130
131
132  if (verbosity.eq.1) then
133     print*,'concoutput_surf 2'
134     CALL SYSTEM_CLOCK(count_clock)
135     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
136  endif
137
138  !*******************************************************************
139  ! Compute air density: sufficiently accurate to take it
140  ! from coarse grid at some time
141  ! Determine center altitude of output layer, and interpolate density
142  ! data to that altitude
143  !*******************************************************************
144
145  mind=memind(2)
146  do kz=1,numzgrid
147    if (kz.eq.1) then
148      halfheight=outheight(1)/2.
149    else
150      halfheight=(outheight(kz)+outheight(kz-1))/2.
151    endif
152    do kzz=2,nz
153      if ((height(kzz-1).lt.halfheight).and. &
154           (height(kzz).gt.halfheight)) goto 46
155    end do
15646   kzz=max(min(kzz,nz),2)
157    dz1=halfheight-height(kzz-1)
158    dz2=height(kzz)-halfheight
159    dz=dz1+dz2
160    do jy=0,numygrid-1
161      do ix=0,numxgrid-1
162        xl=outlon0+real(ix)*dxout
163        yl=outlat0+real(jy)*dyout
164        xl=(xl-xlon0)/dx
165        yl=(yl-ylat0)/dy
166        iix=max(min(nint(xl),nxmin1),0)
167        jjy=max(min(nint(yl),nymin1),0)
168        ! densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
169        !      rho(iix,jjy,kzz-1,2)*dz2)/dz
170        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,mind)*dz1+ &
171             rho(iix,jjy,kzz-1,mind)*dz2)/dz
172      end do
173    end do
174  end do
175
176  do i=1,numreceptor
177    xl=xreceptor(i)
178    yl=yreceptor(i)
179    iix=max(min(nint(xl),nxmin1),0)
180    jjy=max(min(nint(yl),nymin1),0)
181    !densityoutrecept(i)=rho(iix,jjy,1,2)
182    densityoutrecept(i)=rho(iix,jjy,1,mind)
183  end do
184
185
186  ! Output is different for forward and backward simulations
187    do kz=1,numzgrid
188      do jy=0,numygrid-1
189        do ix=0,numxgrid-1
190          if (ldirect.eq.1) then
191            factor3d(ix,jy,kz)=1.e12/volume(ix,jy,kz)/outnum
192          else
193            factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
194          endif
195        end do
196      end do
197    end do
198
199  !*********************************************************************
200  ! Determine the standard deviation of the mean concentration or mixing
201  ! ratio (uncertainty of the output) and the dry and wet deposition
202  !*********************************************************************
203
204  if (verbosity.eq.1) then
205     print*,'concoutput_surf 3 (sd)'
206     CALL SYSTEM_CLOCK(count_clock)
207     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
208  endif
209  gridtotal=0.
210  gridsigmatotal=0.
211  gridtotalunc=0.
212  wetgridtotal=0.
213  wetgridsigmatotal=0.
214  wetgridtotalunc=0.
215  drygridtotal=0.
216  drygridsigmatotal=0.
217  drygridtotalunc=0.
218
219  do ks=1,nspec
220
221  write(anspec,'(i3.3)') ks
222  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
223    if (ldirect.eq.1) then
224      open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_'//adate// &
225           atime//'_'//anspec,form='unformatted')
226    else
227      open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_'//adate// &
228           atime//'_'//anspec,form='unformatted')
229    endif
230     write(unitoutgrid) itime
231   endif
232
233  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
234   open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_'//adate// &
235        atime//'_'//anspec,form='unformatted')
236
237    write(unitoutgridppt) itime
238  endif
239
240  do kp=1,maxpointspec_act
241  do nage=1,nageclass
242
243    do jy=0,numygrid-1
244      do ix=0,numxgrid-1
245
246  ! WET DEPOSITION
247        if ((WETDEP).and.(ldirect.gt.0)) then
248            do l=1,nclassunc
249              auxgrid(l)=wetgridunc0(ix,jy,ks,kp,l,nage)
250            end do
251            call mean(auxgrid,wetgrid(ix,jy), &
252                 wetgridsigma(ix,jy),nclassunc)
253  ! Multiply by number of classes to get total concentration
254            wetgrid(ix,jy)=wetgrid(ix,jy) &
255                 *nclassunc
256            wetgridtotal=wetgridtotal+wetgrid(ix,jy)
257  ! Calculate standard deviation of the mean
258            wetgridsigma(ix,jy)= &
259                 wetgridsigma(ix,jy)* &
260                 sqrt(real(nclassunc))
261            wetgridsigmatotal=wetgridsigmatotal+ &
262                 wetgridsigma(ix,jy)
263        endif
264
265  ! DRY DEPOSITION
266        if ((DRYDEP).and.(ldirect.gt.0)) then
267            do l=1,nclassunc
268              auxgrid(l)=drygridunc0(ix,jy,ks,kp,l,nage)
269            end do
270            call mean(auxgrid,drygrid(ix,jy), &
271                 drygridsigma(ix,jy),nclassunc)
272  ! Multiply by number of classes to get total concentration
273            drygrid(ix,jy)=drygrid(ix,jy)* &
274                 nclassunc
275            drygridtotal=drygridtotal+drygrid(ix,jy)
276  ! Calculate standard deviation of the mean
277            drygridsigma(ix,jy)= &
278                 drygridsigma(ix,jy)* &
279                 sqrt(real(nclassunc))
280125         drygridsigmatotal=drygridsigmatotal+ &
281                 drygridsigma(ix,jy)
282        endif
283
284  ! CONCENTRATION OR MIXING RATIO
285        do kz=1,numzgrid
286            do l=1,nclassunc
287              auxgrid(l)=gridunc(ix,jy,kz,ks,kp,l,nage)
288            end do
289            call mean(auxgrid,grid(ix,jy,kz), &
290                 gridsigma(ix,jy,kz),nclassunc)
291  ! Multiply by number of classes to get total concentration
292            grid(ix,jy,kz)= &
293                 grid(ix,jy,kz)*nclassunc
294            gridtotal=gridtotal+grid(ix,jy,kz)
295  ! Calculate standard deviation of the mean
296            gridsigma(ix,jy,kz)= &
297                 gridsigma(ix,jy,kz)* &
298                 sqrt(real(nclassunc))
299            gridsigmatotal=gridsigmatotal+ &
300                 gridsigma(ix,jy,kz)
301        end do
302      end do
303    end do
304
305
306  !*******************************************************************
307  ! Generate output: may be in concentration (ng/m3) or in mixing
308  ! ratio (ppt) or both
309  ! Output the position and the values alternated multiplied by
310  ! 1 or -1, first line is number of values, number of positions
311  ! For backward simulations, the unit is seconds, stored in grid_time
312  !*******************************************************************
313
314  if (verbosity.eq.1) then
315     print*,'concoutput_surf 4 (output)'
316     CALL SYSTEM_CLOCK(count_clock)
317     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
318  endif
319
320  ! Concentration output
321  !*********************
322
323  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
324
325  if (verbosity.eq.1) then
326     print*,'concoutput_surf (Wet deposition)'
327     CALL SYSTEM_CLOCK(count_clock)
328     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
329  endif
330
331  ! Wet deposition
332         sp_count_i=0
333         sp_count_r=0
334         sp_fact=-1.
335         sp_zer=.true.
336         if ((ldirect.eq.1).and.(WETDEP)) then
337         do jy=0,numygrid-1
338            do ix=0,numxgrid-1
339  ! concentraion greater zero
340              if (wetgrid(ix,jy).gt.smallnum) then
341                 if (sp_zer.eqv..true.) then ! first non zero value
342                    sp_count_i=sp_count_i+1
343                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
344                    sp_zer=.false.
345                    sp_fact=sp_fact*(-1.)
346                 endif
347                 sp_count_r=sp_count_r+1
348                 sparse_dump_r(sp_count_r)= &
349                      sp_fact*1.e12*wetgrid(ix,jy)/area(ix,jy)
350                 sparse_dump_u(sp_count_r)= &
351                      1.e12*wetgridsigma(ix,jy)/area(ix,jy)
352              else ! concentration is zero
353                  sp_zer=.true.
354              endif
355            end do
356         end do
357         else
358            sp_count_i=0
359            sp_count_r=0
360         endif
361         write(unitoutgrid) sp_count_i
362         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
363         write(unitoutgrid) sp_count_r
364         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
365!         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
366
367  if (verbosity.eq.1) then
368     print*,'concoutput_surf (Dry deposition)'
369     CALL SYSTEM_CLOCK(count_clock)
370     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
371  endif
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)/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(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
407
408  if (verbosity.eq.1) then
409     print*,'concoutput_surf (Concentrations)'
410     CALL SYSTEM_CLOCK(count_clock)
411     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
412  endif
413
414  ! Concentrations
415
416  ! surf_only write only 1st layer
417
418         sp_count_i=0
419         sp_count_r=0
420         sp_fact=-1.
421         sp_zer=.true.
422          do kz=1,1
423            do jy=0,numygrid-1
424              do ix=0,numxgrid-1
425                if (grid(ix,jy,kz).gt.smallnum) then
426                  if (sp_zer.eqv..true.) then ! first non zero value
427                    sp_count_i=sp_count_i+1
428                    sparse_dump_i(sp_count_i)= &
429                         ix+jy*numxgrid+kz*numxgrid*numygrid
430                    sp_zer=.false.
431                    sp_fact=sp_fact*(-1.)
432                  endif
433                   sp_count_r=sp_count_r+1
434                   sparse_dump_r(sp_count_r)= &
435                        sp_fact* &
436                        grid(ix,jy,kz)* &
437                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
438  !                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
439  !    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
440                   sparse_dump_u(sp_count_r)= &
441                        gridsigma(ix,jy,kz)* &
442                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
443                else ! concentration is zero
444                  sp_zer=.true.
445                endif
446             end do
447           end do
448         end do
449         write(unitoutgrid) sp_count_i
450         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
451         write(unitoutgrid) sp_count_r
452         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
453!         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
454
455  endif !  concentration output
456
457  ! Mixing ratio output
458  !********************
459
460  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
461
462  ! Wet deposition
463         sp_count_i=0
464         sp_count_r=0
465         sp_fact=-1.
466         sp_zer=.true.
467         if ((ldirect.eq.1).and.(WETDEP)) then
468          do jy=0,numygrid-1
469            do ix=0,numxgrid-1
470                if (wetgrid(ix,jy).gt.smallnum) then
471                  if (sp_zer.eqv..true.) then ! first non zero value
472                    sp_count_i=sp_count_i+1
473                    sparse_dump_i(sp_count_i)= &
474                         ix+jy*numxgrid
475                    sp_zer=.false.
476                    sp_fact=sp_fact*(-1.)
477                 endif
478                 sp_count_r=sp_count_r+1
479                 sparse_dump_r(sp_count_r)= &
480                      sp_fact* &
481                      1.e12*wetgrid(ix,jy)/area(ix,jy)
482                 sparse_dump_u(sp_count_r)= &
483                      1.e12*wetgridsigma(ix,jy)/area(ix,jy)
484              else ! concentration is zero
485                  sp_zer=.true.
486              endif
487            end do
488          end do
489         else
490           sp_count_i=0
491           sp_count_r=0
492         endif
493         write(unitoutgridppt) sp_count_i
494         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
495         write(unitoutgridppt) sp_count_r
496         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
497!         write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
498
499
500  ! Dry deposition
501         sp_count_i=0
502         sp_count_r=0
503         sp_fact=-1.
504         sp_zer=.true.
505         if ((ldirect.eq.1).and.(DRYDEP)) then
506          do jy=0,numygrid-1
507            do ix=0,numxgrid-1
508                if (drygrid(ix,jy).gt.smallnum) then
509                  if (sp_zer.eqv..true.) then ! first non zero value
510                    sp_count_i=sp_count_i+1
511                    sparse_dump_i(sp_count_i)= &
512                         ix+jy*numxgrid
513                    sp_zer=.false.
514                    sp_fact=sp_fact*(-1)
515                 endif
516                 sp_count_r=sp_count_r+1
517                 sparse_dump_r(sp_count_r)= &
518                      sp_fact* &
519                      1.e12*drygrid(ix,jy)/area(ix,jy)
520                 sparse_dump_u(sp_count_r)= &
521                      1.e12*drygridsigma(ix,jy)/area(ix,jy)
522              else ! concentration is zero
523                  sp_zer=.true.
524              endif
525            end do
526          end do
527         else
528           sp_count_i=0
529           sp_count_r=0
530         endif
531         write(unitoutgridppt) sp_count_i
532         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
533         write(unitoutgridppt) sp_count_r
534         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
535!         write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
536
537
538  ! Mixing ratios
539
540  ! surf_only write only 1st layer
541
542         sp_count_i=0
543         sp_count_r=0
544         sp_fact=-1.
545         sp_zer=.true.
546          do kz=1,1
547            do jy=0,numygrid-1
548              do ix=0,numxgrid-1
549                if (grid(ix,jy,kz).gt.smallnum) then
550                  if (sp_zer.eqv..true.) then ! first non zero value
551                    sp_count_i=sp_count_i+1
552                    sparse_dump_i(sp_count_i)= &
553                         ix+jy*numxgrid+kz*numxgrid*numygrid
554                    sp_zer=.false.
555                    sp_fact=sp_fact*(-1.)
556                 endif
557                 sp_count_r=sp_count_r+1
558                 sparse_dump_r(sp_count_r)= &
559                      sp_fact* &
560                      1.e12*grid(ix,jy,kz) &
561                      /volume(ix,jy,kz)/outnum* &
562                      weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
563                 sparse_dump_u(sp_count_r)= &
564                      1.e12*gridsigma(ix,jy,kz)/volume(ix,jy,kz)/ &
565                      outnum*weightair/weightmolar(ks)/ &
566                      densityoutgrid(ix,jy,kz)
567              else ! concentration is zero
568                  sp_zer=.true.
569              endif
570              end do
571            end do
572          end do
573         write(unitoutgridppt) sp_count_i
574         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
575         write(unitoutgridppt) sp_count_r
576         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
577!         write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
578
579      endif ! output for ppt
580
581  end do
582  end do
583
584    close(unitoutgridppt)
585    close(unitoutgrid)
586
587  end do
588
589  if (gridtotal.gt.0.) gridtotalunc=gridsigmatotal/gridtotal
590  if (wetgridtotal.gt.0.) wetgridtotalunc=wetgridsigmatotal/ &
591       wetgridtotal
592  if (drygridtotal.gt.0.) drygridtotalunc=drygridsigmatotal/ &
593       drygridtotal
594
595  ! Dump of receptor concentrations
596
597  if (numreceptor.gt.0 .and. (iout.eq.2 .or. iout.eq.3)  ) then
598    write(unitoutreceptppt) itime
599    do ks=1,nspec
600      write(unitoutreceptppt) (1.e12*creceptor(i,ks)/outnum* &
601           weightair/weightmolar(ks)/densityoutrecept(i),i=1,numreceptor)
602    end do
603  endif
604
605  ! Dump of receptor concentrations
606
607  if (numreceptor.gt.0) then
608    write(unitoutrecept) itime
609    do ks=1,nspec
610      write(unitoutrecept) (1.e12*creceptor(i,ks)/outnum, &
611           i=1,numreceptor)
612    end do
613  endif
614
615
616
617  ! Reinitialization of grid
618  !*************************
619
620  do ks=1,nspec
621    do kp=1,maxpointspec_act
622      do i=1,numreceptor
623        creceptor(i,ks)=0.
624      end do
625      do jy=0,numygrid-1
626        do ix=0,numxgrid-1
627          do l=1,nclassunc
628            do nage=1,nageclass
629              do kz=1,numzgrid
630                gridunc(ix,jy,kz,ks,kp,l,nage)=0.
631              end do
632            end do
633          end do
634        end do
635      end do
636    end do
637  end do
638
639  if (mp_measure_time) call mpif_mtime('rootonly',1)
640 
641end subroutine concoutput_surf
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG