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