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