source: flexpart.git/src/concoutput_surf.f90 @ 3481cc1

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

move license from headers to a different file

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