source: flexpart.git/src/concoutput.f90 @ f9ce123

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

Updated wet depo scheme

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