source: flexpart.git/src/concoutput_surf_nest_mpi.f90

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file was 92fab65, checked in by Ignacio Pisso <ip@…>, 4 years ago

add SPDX-License-Identifier to all .f90 files

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