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