source: flexpart.git/src/concoutput_surf_nest_mpi.f90 @ 4d45639

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

Updated wet depo scheme

  • Property mode set to 100644
File size: 24.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_nest(itime,outnum)
23  !                        i     i
24  !*****************************************************************************
25  !                                                                            *
26  !     Output of the concentration grid and the receptor concentrations.      *
27  !                                                                            *
28  !     Author: A. Stohl                                                       *
29  !                                                                            *
30  !     24 May 1995                                                            *
31  !                                                                            *
32  !     13 April 1999, Major update: if output size is smaller, dump output    *
33  !                    in sparse matrix format; additional output of           *
34  !                    uncertainty                                             *
35  !                                                                            *
36  !     05 April 2000, Major update: output of age classes; output for backward*
37  !                    runs is time spent in grid cell times total mass of     *
38  !                    species.                                                *
39  !                                                                            *
40  !     17 February 2002, Appropriate dimensions for backward and forward runs *
41  !                       are now specified in file par_mod                    *
42  !                                                                            *
43  !     June 2006, write grid in sparse matrix with a single write command     *
44  !                in order to save disk space                                 *
45  !                                                                            *
46  !     2008 new sparse matrix format                                          *
47  !                                                                            *
48  !     2014 eso: MPI version. Only called by root process                     *
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  use mpi_mod
66
67  implicit none
68
69  real(kind=dp) :: jul
70  integer :: itime,i,ix,jy,kz,ks,kp,l,iix,jjy,kzz,nage,jjjjmmdd,ihmmss
71  integer :: sp_count_i,sp_count_r
72  real :: sp_fact
73  real :: outnum,densityoutrecept(maxreceptor),xl,yl
74
75  !real densityoutgrid(0:numxgrid-1,0:numygrid-1,numzgrid),
76  !    +grid(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,maxpointspec_act,
77  !    +    maxageclass)
78  !real wetgrid(0:numxgrid-1,0:numygrid-1,maxspec,maxpointspec_act,
79  !    +       maxageclass)
80  !real drygrid(0:numxgrid-1,0:numygrid-1,maxspec,
81  !    +       maxpointspec_act,maxageclass)
82  !real gridsigma(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,
83  !    +       maxpointspec_act,maxageclass),
84  !    +     drygridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
85  !    +     maxpointspec_act,maxageclass),
86  !    +     wetgridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
87  !    +     maxpointspec_act,maxageclass)
88  !real factor(0:numxgrid-1,0:numygrid-1,numzgrid)
89  !real sparse_dump_r(numxgrid*numygrid*numzgrid)
90  !integer sparse_dump_i(numxgrid*numygrid*numzgrid)
91
92  !real sparse_dump_u(numxgrid*numygrid*numzgrid)
93  real :: auxgrid(nclassunc)
94  real :: halfheight,dz,dz1,dz2,tot_mu(maxspec,maxpointspec_act)
95  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
96  real,parameter :: weightair=28.97
97  logical :: sp_zer
98  character :: adate*8,atime*6
99  character(len=3) :: anspec
100  integer :: mind
101! mind        eso:added to get consistent results between 2&3-fields versions
102
103! Measure execution time
104  if (mp_measure_time) then
105    call cpu_time(mp_root_time_beg)
106    mp_root_wtime_beg = mpi_wtime()
107  end if
108
109  if (verbosity.eq.1) then
110     print*,'inside concoutput_surf '
111     CALL SYSTEM_CLOCK(count_clock)
112     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
113  endif
114
115  ! Determine current calendar date, needed for the file name
116  !**********************************************************
117
118  jul=bdate+real(itime,kind=dp)/86400._dp
119  call caldate(jul,jjjjmmdd,ihmmss)
120  write(adate,'(i8.8)') jjjjmmdd
121  write(atime,'(i6.6)') ihmmss
122
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,numygridn-1
168      do ix=0,numxgridn-1
169        xl=outlon0n+real(ix)*dxoutn
170        yl=outlat0n+real(jy)*dyoutn
171        xl=(xl-xlon0)/dx
172        yl=(yl-ylat0)/dy
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,numygridn-1
196        do ix=0,numxgridn-1
197          if (ldirect.eq.1) then
198            factor3d(ix,jy,kz)=1.e12/volumen(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  do ks=1,nspec
212
213  write(anspec,'(i3.3)') ks
214  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
215    if (ldirect.eq.1) then
216      open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_nest_' &
217           //adate// &
218           atime//'_'//anspec,form='unformatted')
219    else
220      open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_nest_' &
221           //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_nest_' &
229        //adate// &
230        atime//'_'//anspec,form='unformatted')
231
232    write(unitoutgridppt) itime
233  endif
234
235  do kp=1,maxpointspec_act
236  do nage=1,nageclass
237
238    do jy=0,numygridn-1
239      do ix=0,numxgridn-1
240
241  ! WET DEPOSITION
242        if ((WETDEP).and.(ldirect.gt.0)) then
243            do l=1,nclassunc
244              auxgrid(l)=wetgriduncn0(ix,jy,ks,kp,l,nage)
245            end do
246            call mean(auxgrid,wetgrid(ix,jy), &
247                 wetgridsigma(ix,jy),nclassunc)
248  ! Multiply by number of classes to get total concentration
249            wetgrid(ix,jy)=wetgrid(ix,jy) &
250                 *nclassunc
251  ! Calculate standard deviation of the mean
252            wetgridsigma(ix,jy)= &
253                 wetgridsigma(ix,jy)* &
254                 sqrt(real(nclassunc))
255        endif
256
257  ! DRY DEPOSITION
258        if ((DRYDEP).and.(ldirect.gt.0)) then
259            do l=1,nclassunc
260              auxgrid(l)=drygriduncn0(ix,jy,ks,kp,l,nage)
261            end do
262            call mean(auxgrid,drygrid(ix,jy), &
263                 drygridsigma(ix,jy),nclassunc)
264  ! Multiply by number of classes to get total concentration
265            drygrid(ix,jy)=drygrid(ix,jy)* &
266                 nclassunc
267  ! Calculate standard deviation of the mean
268            drygridsigma(ix,jy)= &
269                 drygridsigma(ix,jy)* &
270                 sqrt(real(nclassunc))
271        endif
272
273  ! CONCENTRATION OR MIXING RATIO
274        do kz=1,numzgrid
275            do l=1,nclassunc
276              auxgrid(l)=griduncn(ix,jy,kz,ks,kp,l,nage)
277            end do
278            call mean(auxgrid,grid(ix,jy,kz), &
279                 gridsigma(ix,jy,kz),nclassunc)
280  ! Multiply by number of classes to get total concentration
281            grid(ix,jy,kz)= &
282                 grid(ix,jy,kz)*nclassunc
283  ! Calculate standard deviation of the mean
284            gridsigma(ix,jy,kz)= &
285                 gridsigma(ix,jy,kz)* &
286                 sqrt(real(nclassunc))
287        end do
288      end do
289    end do
290
291
292  !*******************************************************************
293  ! Generate output: may be in concentration (ng/m3) or in mixing
294  ! ratio (ppt) or both
295  ! Output the position and the values alternated multiplied by
296  ! 1 or -1, first line is number of values, number of positions
297  ! For backward simulations, the unit is seconds, stored in grid_time
298  !*******************************************************************
299
300  ! Concentration output
301  !*********************
302  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
303
304  ! Wet deposition
305         sp_count_i=0
306         sp_count_r=0
307         sp_fact=-1.
308         sp_zer=.true.
309         if ((ldirect.eq.1).and.(WETDEP)) then
310         do jy=0,numygridn-1
311            do ix=0,numxgridn-1
312  !oncentraion greater zero
313              if (wetgrid(ix,jy).gt.smallnum) then
314                 if (sp_zer.eqv..true.) then ! first non zero value
315                    sp_count_i=sp_count_i+1
316                    sparse_dump_i(sp_count_i)=ix+jy*numxgridn
317                    sp_zer=.false.
318                    sp_fact=sp_fact*(-1.)
319                 endif
320                 sp_count_r=sp_count_r+1
321                 sparse_dump_r(sp_count_r)= &
322                      sp_fact*1.e12*wetgrid(ix,jy)/arean(ix,jy)
323                 sparse_dump_u(sp_count_r)= &
324                      1.e12*wetgridsigma(ix,jy)/area(ix,jy)
325              else ! concentration is zero
326                  sp_zer=.true.
327              endif
328            end do
329         end do
330         else
331            sp_count_i=0
332            sp_count_r=0
333         endif
334         write(unitoutgrid) sp_count_i
335         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
336         write(unitoutgrid) sp_count_r
337         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
338         write(unitoutgrid) sp_count_r
339         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
340
341  ! Dry deposition
342         sp_count_i=0
343         sp_count_r=0
344         sp_fact=-1.
345         sp_zer=.true.
346         if ((ldirect.eq.1).and.(DRYDEP)) then
347          do jy=0,numygridn-1
348            do ix=0,numxgridn-1
349              if (drygrid(ix,jy).gt.smallnum) then
350                 if (sp_zer.eqv..true.) then ! first non zero value
351                    sp_count_i=sp_count_i+1
352                    sparse_dump_i(sp_count_i)=ix+jy*numxgridn
353                    sp_zer=.false.
354                    sp_fact=sp_fact*(-1.)
355                 endif
356                 sp_count_r=sp_count_r+1
357                 sparse_dump_r(sp_count_r)= &
358                      sp_fact* &
359                      1.e12*drygrid(ix,jy)/arean(ix,jy)
360                 sparse_dump_u(sp_count_r)= &
361                      1.e12*drygridsigma(ix,jy)/area(ix,jy)
362              else ! concentration is zero
363                  sp_zer=.true.
364              endif
365            end do
366          end do
367         else
368            sp_count_i=0
369            sp_count_r=0
370         endif
371         write(unitoutgrid) sp_count_i
372         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
373         write(unitoutgrid) sp_count_r
374         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
375         write(unitoutgrid) sp_count_r
376         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
377
378
379
380  ! Concentrations
381
382  ! if surf_only write only 1st layer
383
384         if(surf_only.eq.1) then
385         sp_count_i=0
386         sp_count_r=0
387         sp_fact=-1.
388         sp_zer=.true.
389          do kz=1,1
390            do jy=0,numygridn-1
391              do ix=0,numxgridn-1
392                if (grid(ix,jy,kz).gt.smallnum) then
393                  if (sp_zer.eqv..true.) then ! first non zero value
394                    sp_count_i=sp_count_i+1
395                    sparse_dump_i(sp_count_i)= &
396                         ix+jy*numxgridn+kz*numxgridn*numygridn
397                    sp_zer=.false.
398                    sp_fact=sp_fact*(-1.)
399                   endif
400                   sp_count_r=sp_count_r+1
401                   sparse_dump_r(sp_count_r)= &
402                        sp_fact* &
403                        grid(ix,jy,kz)* &
404                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
405  !                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
406  !    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
407                   sparse_dump_u(sp_count_r)= &
408                        gridsigma(ix,jy,kz)* &
409                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
410              else ! concentration is zero
411                  sp_zer=.true.
412              endif
413              end do
414            end do
415          end do
416         write(unitoutgrid) sp_count_i
417         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
418         write(unitoutgrid) sp_count_r
419         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
420         write(unitoutgrid) sp_count_r
421         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
422         else
423
424  ! write full vertical resolution
425
426         sp_count_i=0
427         sp_count_r=0
428         sp_fact=-1.
429         sp_zer=.true.
430          do kz=1,numzgrid
431            do jy=0,numygridn-1
432              do ix=0,numxgridn-1
433                if (grid(ix,jy,kz).gt.smallnum) then
434                  if (sp_zer.eqv..true.) then ! first non zero value
435                    sp_count_i=sp_count_i+1
436                    sparse_dump_i(sp_count_i)= &
437                         ix+jy*numxgridn+kz*numxgridn*numygridn
438                    sp_zer=.false.
439                    sp_fact=sp_fact*(-1.)
440                   endif
441                   sp_count_r=sp_count_r+1
442                   sparse_dump_r(sp_count_r)= &
443                        sp_fact* &
444                        grid(ix,jy,kz)* &
445                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
446  !                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
447  !    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
448                   sparse_dump_u(sp_count_r)= &
449                        gridsigma(ix,jy,kz)* &
450                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
451              else ! concentration is zero
452                  sp_zer=.true.
453              endif
454              end do
455            end do
456          end do
457         write(unitoutgrid) sp_count_i
458         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
459         write(unitoutgrid) sp_count_r
460         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
461         write(unitoutgrid) sp_count_r
462         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
463         endif ! surf_only
464
465
466    endif !  concentration output
467
468  ! Mixing ratio output
469  !********************
470
471  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
472
473  ! Wet deposition
474         sp_count_i=0
475         sp_count_r=0
476         sp_fact=-1.
477         sp_zer=.true.
478         if ((ldirect.eq.1).and.(WETDEP)) then
479          do jy=0,numygridn-1
480            do ix=0,numxgridn-1
481                if (wetgrid(ix,jy).gt.smallnum) then
482                  if (sp_zer.eqv..true.) then ! first non zero value
483                    sp_count_i=sp_count_i+1
484                    sparse_dump_i(sp_count_i)= &
485                         ix+jy*numxgridn
486                    sp_zer=.false.
487                    sp_fact=sp_fact*(-1.)
488                 endif
489                 sp_count_r=sp_count_r+1
490                 sparse_dump_r(sp_count_r)= &
491                      sp_fact* &
492                      1.e12*wetgrid(ix,jy)/arean(ix,jy)
493                 sparse_dump_u(sp_count_r)= &
494                      1.e12*wetgridsigma(ix,jy)/area(ix,jy)
495              else ! concentration is zero
496                  sp_zer=.true.
497              endif
498            end do
499          end do
500         else
501           sp_count_i=0
502           sp_count_r=0
503         endif
504         write(unitoutgridppt) sp_count_i
505         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
506         write(unitoutgridppt) sp_count_r
507         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
508         write(unitoutgridppt) sp_count_r
509         write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
510
511
512  ! Dry deposition
513         sp_count_i=0
514         sp_count_r=0
515         sp_fact=-1.
516         sp_zer=.true.
517         if ((ldirect.eq.1).and.(DRYDEP)) then
518          do jy=0,numygridn-1
519            do ix=0,numxgridn-1
520                if (drygrid(ix,jy).gt.smallnum) then
521                  if (sp_zer.eqv..true.) then ! first non zero value
522                    sp_count_i=sp_count_i+1
523                    sparse_dump_i(sp_count_i)= &
524                         ix+jy*numxgridn
525                    sp_zer=.false.
526                    sp_fact=sp_fact*(-1)
527                 endif
528                 sp_count_r=sp_count_r+1
529                 sparse_dump_r(sp_count_r)= &
530                      sp_fact* &
531                      1.e12*drygrid(ix,jy)/arean(ix,jy)
532                 sparse_dump_u(sp_count_r)= &
533                      1.e12*drygridsigma(ix,jy)/area(ix,jy)
534              else ! concentration is zero
535                  sp_zer=.true.
536              endif
537            end do
538          end do
539         else
540           sp_count_i=0
541           sp_count_r=0
542         endif
543         write(unitoutgridppt) sp_count_i
544         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
545         write(unitoutgridppt) sp_count_r
546         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
547         write(unitoutgridppt) sp_count_r
548         write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
549
550
551  ! Mixing ratios
552
553    ! if surf_only write only 1st layer
554
555         if(surf_only.eq.1) then
556         sp_count_i=0
557         sp_count_r=0
558         sp_fact=-1.
559         sp_zer=.true.
560          do kz=1,1
561            do jy=0,numygridn-1
562              do ix=0,numxgridn-1
563                if (grid(ix,jy,kz).gt.smallnum) then
564                  if (sp_zer.eqv..true.) then ! first non zero value
565                    sp_count_i=sp_count_i+1
566                    sparse_dump_i(sp_count_i)= &
567                         ix+jy*numxgridn+kz*numxgridn*numygridn
568                    sp_zer=.false.
569                    sp_fact=sp_fact*(-1.)
570                 endif
571                 sp_count_r=sp_count_r+1
572                 sparse_dump_r(sp_count_r)= &
573                      sp_fact* &
574                      1.e12*grid(ix,jy,kz) &
575                      /volumen(ix,jy,kz)/outnum* &
576                      weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
577                 sparse_dump_u(sp_count_r)= &
578                      1.e12*gridsigma(ix,jy,kz)/volumen(ix,jy,kz)/ &
579                      outnum*weightair/weightmolar(ks)/ &
580                      densityoutgrid(ix,jy,kz)
581              else ! concentration is zero
582                  sp_zer=.true.
583              endif
584              end do
585            end do
586          end do
587         write(unitoutgridppt) sp_count_i
588         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
589         write(unitoutgridppt) sp_count_r
590         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
591         write(unitoutgridppt) sp_count_r
592         write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
593         else
594
595  ! write full vertical resolution
596
597         sp_count_i=0
598         sp_count_r=0
599         sp_fact=-1.
600         sp_zer=.true.
601          do kz=1,numzgrid
602            do jy=0,numygridn-1
603              do ix=0,numxgridn-1
604                if (grid(ix,jy,kz).gt.smallnum) then
605                  if (sp_zer.eqv..true.) then ! first non zero value
606                    sp_count_i=sp_count_i+1
607                    sparse_dump_i(sp_count_i)= &
608                         ix+jy*numxgridn+kz*numxgridn*numygridn
609                    sp_zer=.false.
610                    sp_fact=sp_fact*(-1.)
611                 endif
612                 sp_count_r=sp_count_r+1
613                 sparse_dump_r(sp_count_r)= &
614                      sp_fact* &
615                      1.e12*grid(ix,jy,kz) &
616                      /volumen(ix,jy,kz)/outnum* &
617                      weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
618                 sparse_dump_u(sp_count_r)= &
619                      1.e12*gridsigma(ix,jy,kz)/volumen(ix,jy,kz)/ &
620                      outnum*weightair/weightmolar(ks)/ &
621                      densityoutgrid(ix,jy,kz)
622              else ! concentration is zero
623                  sp_zer=.true.
624              endif
625              end do
626            end do
627          end do
628         write(unitoutgridppt) sp_count_i
629         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
630         write(unitoutgridppt) sp_count_r
631         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
632         write(unitoutgridppt) sp_count_r
633         write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
634         endif ! surf_only
635
636      endif ! output for ppt
637
638  end do
639  end do
640
641    close(unitoutgridppt)
642    close(unitoutgrid)
643
644  end do
645
646
647
648  ! Reinitialization of grid
649  !*************************
650
651  do ks=1,nspec
652  do kp=1,maxpointspec_act
653    do i=1,numreceptor
654      creceptor(i,ks)=0.
655    end do
656    do jy=0,numygridn-1
657      do ix=0,numxgridn-1
658        do l=1,nclassunc
659          do nage=1,nageclass
660            do kz=1,numzgrid
661              griduncn(ix,jy,kz,ks,kp,l,nage)=0.
662            end do
663          end do
664        end do
665      end do
666    end do
667  end do
668  end do
669
670  if (mp_measure_time) then
671    call cpu_time(mp_root_time_end)
672    mp_root_wtime_end = mpi_wtime()
673    mp_root_time_total = mp_root_time_total + (mp_root_time_end - mp_root_time_beg)
674    mp_root_wtime_total = mp_root_wtime_total + (mp_root_wtime_end - mp_root_wtime_beg)
675  end if
676
677end subroutine concoutput_surf_nest
678
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG