source: flexpart.git/src/concoutput.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: 25.3 KB
RevLine 
[92fab65]1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
[332fbbd]3
[e200b7a]4subroutine concoutput(itime,outnum,gridtotalunc,wetgridtotalunc, &
[b069789]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!*****************************************************************************
[e200b7a]41
42  use unc_mod
43  use point_mod
44  use outg_mod
45  use par_mod
46  use com_mod
[6a678e3]47  use mean_mod
[e200b7a]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
[2eefa58]56! RLT
57  real :: densitydryrecept(maxreceptor)
58  real :: factor_dryrecept(maxreceptor)
[e200b7a]59
[b069789]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)
[6a678e3]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
[e200b7a]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
[b069789]86  logical,save :: init=.true.
[e200b7a]87  character :: adate*8,atime*6
88  character(len=3) :: anspec
[5f9d14a]89  integer :: mind
90! mind        eso:added to ensure identical results between 2&3-fields versions
[b069789]91  character(LEN=8),save :: file_stat='REPLACE'
92  logical :: ldates_file
[2eefa58]93  logical :: lexist
[b069789]94  integer :: ierr
95  character(LEN=100) :: dates_char
[e200b7a]96
[b069789]97! Determine current calendar date, needed for the file name
98!**********************************************************
[e200b7a]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
[b069789]104
105! Overwrite existing dates file on first call, later append to it
106! This fixes a bug where the dates file kept growing across multiple runs
107
[6b22af9]108! If 'dates' file exists in output directory, make a backup
[b069789]109  inquire(file=path(2)(1:length(2))//'dates', exist=ldates_file)
110  if (ldates_file.and.init) then
111    open(unit=unitdates, file=path(2)(1:length(2))//'dates',form='formatted', &
112         &access='sequential', status='old', action='read', iostat=ierr)
113    open(unit=unittmp, file=path(2)(1:length(2))//'dates.bak', access='sequential', &
114         &status='replace', action='write', form='formatted', iostat=ierr)
115    do while (.true.)
116      read(unitdates, '(a)', iostat=ierr) dates_char
117      if (ierr.ne.0) exit
118      write(unit=unittmp, fmt='(a)', iostat=ierr, advance='yes') trim(dates_char)
119    end do
120    close(unit=unitdates)
121    close(unit=unittmp)
122  end if
123
124  open(unitdates,file=path(2)(1:length(2))//'dates', ACCESS='APPEND', STATUS=file_stat)
[e200b7a]125  write(unitdates,'(a)') adate//atime
[f13406c]126  close(unitdates) 
[e200b7a]127
[20963b1]128  ! Overwrite existing dates file on first call, later append to it
129  ! This fixes a bug where the dates file kept growing across multiple runs
[478e9e6]130  IF (init) THEN
131    file_stat='OLD'
132    init=.false.
133  END IF
134
135
[b069789]136! For forward simulations, output fields have dimension MAXSPEC,
137! for backward simulations, output fields have dimension MAXPOINT.
138! Thus, make loops either about nspec, or about numpoint
139!*****************************************************************
[e200b7a]140
141
142  if (ldirect.eq.1) then
143    do ks=1,nspec
144      do kp=1,maxpointspec_act
145        tot_mu(ks,kp)=1
146      end do
147    end do
148  else
149    do ks=1,nspec
150      do kp=1,maxpointspec_act
151        tot_mu(ks,kp)=xmass(kp,ks)
152      end do
153    end do
154  endif
155
156
[b069789]157!*******************************************************************
158! Compute air density: sufficiently accurate to take it
159! from coarse grid at some time
160! Determine center altitude of output layer, and interpolate density
161! data to that altitude
162!*******************************************************************
[e200b7a]163
[5f9d14a]164  mind=memind(2)
[e200b7a]165  do kz=1,numzgrid
166    if (kz.eq.1) then
167      halfheight=outheight(1)/2.
168    else
169      halfheight=(outheight(kz)+outheight(kz-1))/2.
170    endif
171    do kzz=2,nz
172      if ((height(kzz-1).lt.halfheight).and. &
173           (height(kzz).gt.halfheight)) goto 46
174    end do
[b069789]17546  kzz=max(min(kzz,nz),2)
[e200b7a]176    dz1=halfheight-height(kzz-1)
177    dz2=height(kzz)-halfheight
178    dz=dz1+dz2
179    do jy=0,numygrid-1
180      do ix=0,numxgrid-1
181        xl=outlon0+real(ix)*dxout
182        yl=outlat0+real(jy)*dyout
183        xl=(xl-xlon0)/dx
[f13406c]184        yl=(yl-ylat0)/dy !v9.1.1
[e200b7a]185        iix=max(min(nint(xl),nxmin1),0)
186        jjy=max(min(nint(yl),nymin1),0)
[b069789]187! densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
188!      rho(iix,jjy,kzz-1,2)*dz2)/dz
[5f9d14a]189        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,mind)*dz1+ &
190             rho(iix,jjy,kzz-1,mind)*dz2)/dz
[2eefa58]191! RLT
192        densitydrygrid(ix,jy,kz)=(rho_dry(iix,jjy,kzz,mind)*dz1+ &
193             rho_dry(iix,jjy,kzz-1,mind)*dz2)/dz 
[e200b7a]194      end do
195    end do
196  end do
197
[5f9d14a]198  do i=1,numreceptor
199    xl=xreceptor(i)
200    yl=yreceptor(i)
201    iix=max(min(nint(xl),nxmin1),0)
202    jjy=max(min(nint(yl),nymin1),0)
[b069789]203!densityoutrecept(i)=rho(iix,jjy,1,2)
[5f9d14a]204    densityoutrecept(i)=rho(iix,jjy,1,mind)
[2eefa58]205! RLT
206    densitydryrecept(i)=rho_dry(iix,jjy,1,mind)
[5f9d14a]207  end do
[e200b7a]208
[2eefa58]209! RLT
210! conversion factor for output relative to dry air
211  factor_drygrid=densityoutgrid/densitydrygrid
212  factor_dryrecept=densityoutrecept/densitydryrecept
[e200b7a]213
[b069789]214! Output is different for forward and backward simulations
215  do kz=1,numzgrid
216    do jy=0,numygrid-1
217      do ix=0,numxgrid-1
218        if (ldirect.eq.1) then
219          factor3d(ix,jy,kz)=1.e12/volume(ix,jy,kz)/outnum
220        else
221          factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
222        endif
[e200b7a]223      end do
224    end do
[b069789]225  end do
[e200b7a]226
[b069789]227!*********************************************************************
228! Determine the standard deviation of the mean concentration or mixing
229! ratio (uncertainty of the output) and the dry and wet deposition
230!*********************************************************************
[e200b7a]231
232  gridtotal=0.
233  gridsigmatotal=0.
234  gridtotalunc=0.
235  wetgridtotal=0.
236  wetgridsigmatotal=0.
237  wetgridtotalunc=0.
238  drygridtotal=0.
239  drygridsigmatotal=0.
240  drygridtotalunc=0.
241
242  do ks=1,nspec
243
[b069789]244    write(anspec,'(i3.3)') ks
[e200b7a]245
[54cbd6c]246    if (DRYBKDEP.or.WETBKDEP) then !scavdep output
247      if (DRYBKDEP) &
248      open(unitoutgrid,file=path(2)(1:length(2))//'grid_drydep_'//adate// &
[b069789]249           atime//'_'//anspec,form='unformatted')
[54cbd6c]250      if (WETBKDEP) &
251      open(unitoutgrid,file=path(2)(1:length(2))//'grid_wetdep_'//adate// &
[462f74b]252           atime//'_'//anspec,form='unformatted')
253      write(unitoutgrid) itime
[54cbd6c]254    else
255      if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
256        if (ldirect.eq.1) then
257          open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_'//adate// &
[b069789]258             atime//'_'//anspec,form='unformatted')
[54cbd6c]259        else
260          open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_'//adate// &
[b069789]261             atime//'_'//anspec,form='unformatted')
[54cbd6c]262        endif
263        write(unitoutgrid) itime
[b069789]264      endif
[54cbd6c]265      if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
266        open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_'//adate// &
[b069789]267           atime//'_'//anspec,form='unformatted')
[54cbd6c]268        write(unitoutgridppt) itime
269      endif
270    endif ! if deposition output
[e200b7a]271
[b069789]272    do kp=1,maxpointspec_act
273      do nage=1,nageclass
[e200b7a]274
[b069789]275        do jy=0,numygrid-1
276          do ix=0,numxgrid-1
[e200b7a]277
[b069789]278! WET DEPOSITION
279            if ((WETDEP).and.(ldirect.gt.0)) then
280              do l=1,nclassunc
281                auxgrid(l)=wetgridunc(ix,jy,ks,kp,l,nage)
282              end do
283              call mean(auxgrid,wetgrid(ix,jy), &
284                   wetgridsigma(ix,jy),nclassunc)
285! Multiply by number of classes to get total concentration
286              wetgrid(ix,jy)=wetgrid(ix,jy) &
287                   *nclassunc
288              wetgridtotal=wetgridtotal+wetgrid(ix,jy)
289! Calculate standard deviation of the mean
290              wetgridsigma(ix,jy)= &
291                   wetgridsigma(ix,jy)* &
292                   sqrt(real(nclassunc))
293              wetgridsigmatotal=wetgridsigmatotal+ &
294                   wetgridsigma(ix,jy)
295            endif
296
297! DRY DEPOSITION
298            if ((DRYDEP).and.(ldirect.gt.0)) then
299              do l=1,nclassunc
300                auxgrid(l)=drygridunc(ix,jy,ks,kp,l,nage)
301              end do
302              call mean(auxgrid,drygrid(ix,jy), &
303                   drygridsigma(ix,jy),nclassunc)
304! Multiply by number of classes to get total concentration
305              drygrid(ix,jy)=drygrid(ix,jy)* &
306                   nclassunc
307              drygridtotal=drygridtotal+drygrid(ix,jy)
308! Calculate standard deviation of the mean
309              drygridsigma(ix,jy)= &
310                   drygridsigma(ix,jy)* &
311                   sqrt(real(nclassunc))
[20963b1]312              drygridsigmatotal=drygridsigmatotal+ &
[b069789]313                   drygridsigma(ix,jy)
314            endif
315
316! CONCENTRATION OR MIXING RATIO
317            do kz=1,numzgrid
318              do l=1,nclassunc
319                auxgrid(l)=gridunc(ix,jy,kz,ks,kp,l,nage)
320              end do
321              call mean(auxgrid,grid(ix,jy,kz), &
322                   gridsigma(ix,jy,kz),nclassunc)
323! Multiply by number of classes to get total concentration
324              grid(ix,jy,kz)= &
325                   grid(ix,jy,kz)*nclassunc
326              gridtotal=gridtotal+grid(ix,jy,kz)
327! Calculate standard deviation of the mean
328              gridsigma(ix,jy,kz)= &
329                   gridsigma(ix,jy,kz)* &
330                   sqrt(real(nclassunc))
331              gridsigmatotal=gridsigmatotal+ &
332                   gridsigma(ix,jy,kz)
[e200b7a]333            end do
[b069789]334          end do
[e200b7a]335        end do
336
337
338
339
[b069789]340!*******************************************************************
341! Generate output: may be in concentration (ng/m3) or in mixing
342! ratio (ppt) or both
343! Output the position and the values alternated multiplied by
344! 1 or -1, first line is number of values, number of positions
345! For backward simulations, the unit is seconds, stored in grid_time
346!*******************************************************************
[e200b7a]347
[b069789]348! Concentration output
349!*********************
[2eefa58]350        if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
[e200b7a]351
[b069789]352! Wet deposition
353          sp_count_i=0
354          sp_count_r=0
355          sp_fact=-1.
356          sp_zer=.true.
357          if ((ldirect.eq.1).and.(WETDEP)) then
358            do jy=0,numygrid-1
359              do ix=0,numxgrid-1
360!oncentraion greater zero
361                if (wetgrid(ix,jy).gt.smallnum) then
362                  if (sp_zer.eqv..true.) then ! first non zero value
[e200b7a]363                    sp_count_i=sp_count_i+1
364                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
365                    sp_zer=.false.
366                    sp_fact=sp_fact*(-1.)
[b069789]367                  endif
368                  sp_count_r=sp_count_r+1
369                  sparse_dump_r(sp_count_r)= &
370                       sp_fact*1.e12*wetgrid(ix,jy)/area(ix,jy)
371!                sparse_dump_u(sp_count_r)=
372!+                1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
373                else ! concentration is zero
[e200b7a]374                  sp_zer=.true.
[b069789]375                endif
376              end do
[e200b7a]377            end do
[b069789]378          else
[e200b7a]379            sp_count_i=0
380            sp_count_r=0
[b069789]381          endif
382          write(unitoutgrid) sp_count_i
383          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
384          write(unitoutgrid) sp_count_r
385          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
386!       write(unitoutgrid) sp_count_u
387!       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
388
389! Dry deposition
390          sp_count_i=0
391          sp_count_r=0
392          sp_fact=-1.
393          sp_zer=.true.
394          if ((ldirect.eq.1).and.(DRYDEP)) then
395            do jy=0,numygrid-1
396              do ix=0,numxgrid-1
397                if (drygrid(ix,jy).gt.smallnum) then
398                  if (sp_zer.eqv..true.) then ! first non zero value
[e200b7a]399                    sp_count_i=sp_count_i+1
400                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
401                    sp_zer=.false.
402                    sp_fact=sp_fact*(-1.)
[b069789]403                  endif
404                  sp_count_r=sp_count_r+1
405                  sparse_dump_r(sp_count_r)= &
406                       sp_fact* &
407                       1.e12*drygrid(ix,jy)/area(ix,jy)
408!                sparse_dump_u(sp_count_r)=
409!+                1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
410                else ! concentration is zero
[e200b7a]411                  sp_zer=.true.
[b069789]412                endif
413              end do
[e200b7a]414            end do
[b069789]415          else
[e200b7a]416            sp_count_i=0
417            sp_count_r=0
[b069789]418          endif
419          write(unitoutgrid) sp_count_i
420          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
421          write(unitoutgrid) sp_count_r
422          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
423!       write(*,*) sp_count_u
424!       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
[e200b7a]425
426
427
[b069789]428! Concentrations
429          sp_count_i=0
430          sp_count_r=0
431          sp_fact=-1.
432          sp_zer=.true.
[e200b7a]433          do kz=1,numzgrid
434            do jy=0,numygrid-1
435              do ix=0,numxgrid-1
436                if (grid(ix,jy,kz).gt.smallnum) then
437                  if (sp_zer.eqv..true.) then ! first non zero value
438                    sp_count_i=sp_count_i+1
439                    sparse_dump_i(sp_count_i)= &
440                         ix+jy*numxgrid+kz*numxgrid*numygrid
441                    sp_zer=.false.
442                    sp_fact=sp_fact*(-1.)
[b069789]443                  endif
444                  sp_count_r=sp_count_r+1
[08a38b5]445                  if (lparticlecountoutput) then
446                    sparse_dump_r(sp_count_r)= &
447                         sp_fact* &
448                         grid(ix,jy,kz)
449                  else
450                    sparse_dump_r(sp_count_r)= &
451                         sp_fact* &
452                         grid(ix,jy,kz)* &
453                         factor3d(ix,jy,kz)/tot_mu(ks,kp)
454                  end if
455
456
[b069789]457!                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
458!    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
459!                sparse_dump_u(sp_count_r)=
460!+               ,gridsigma(ix,jy,kz,ks,kp,nage)*
461!+               factor(ix,jy,kz)/tot_mu(ks,kp)
462                else ! concentration is zero
[e200b7a]463                  sp_zer=.true.
[b069789]464                endif
[e200b7a]465              end do
466            end do
467          end do
[b069789]468          write(unitoutgrid) sp_count_i
469          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
470          write(unitoutgrid) sp_count_r
471          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
472!       write(unitoutgrid) sp_count_u
473!       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
[e200b7a]474
475
476
[b069789]477        endif !  concentration output
[e200b7a]478
[b069789]479! Mixing ratio output
480!********************
[e200b7a]481
[b069789]482        if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
[e200b7a]483
[b069789]484! Wet deposition
485          sp_count_i=0
486          sp_count_r=0
487          sp_fact=-1.
488          sp_zer=.true.
489          if ((ldirect.eq.1).and.(WETDEP)) then
490            do jy=0,numygrid-1
491              do ix=0,numxgrid-1
[e200b7a]492                if (wetgrid(ix,jy).gt.smallnum) then
493                  if (sp_zer.eqv..true.) then ! first non zero value
494                    sp_count_i=sp_count_i+1
495                    sparse_dump_i(sp_count_i)= &
496                         ix+jy*numxgrid
497                    sp_zer=.false.
498                    sp_fact=sp_fact*(-1.)
[b069789]499                  endif
500                  sp_count_r=sp_count_r+1
501                  sparse_dump_r(sp_count_r)= &
502                       sp_fact* &
503                       1.e12*wetgrid(ix,jy)/area(ix,jy)
504!                sparse_dump_u(sp_count_r)=
505!    +            ,1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
506                else ! concentration is zero
[e200b7a]507                  sp_zer=.true.
[b069789]508                endif
509              end do
[e200b7a]510            end do
[b069789]511          else
512            sp_count_i=0
513            sp_count_r=0
514          endif
515          write(unitoutgridppt) sp_count_i
516          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
517          write(unitoutgridppt) sp_count_r
518          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
519!       write(unitoutgridppt) sp_count_u
520!       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
521
522
523! Dry deposition
524          sp_count_i=0
525          sp_count_r=0
526          sp_fact=-1.
527          sp_zer=.true.
528          if ((ldirect.eq.1).and.(DRYDEP)) then
529            do jy=0,numygrid-1
530              do ix=0,numxgrid-1
[e200b7a]531                if (drygrid(ix,jy).gt.smallnum) then
532                  if (sp_zer.eqv..true.) then ! first non zero value
533                    sp_count_i=sp_count_i+1
534                    sparse_dump_i(sp_count_i)= &
535                         ix+jy*numxgrid
536                    sp_zer=.false.
537                    sp_fact=sp_fact*(-1)
[b069789]538                  endif
539                  sp_count_r=sp_count_r+1
540                  sparse_dump_r(sp_count_r)= &
541                       sp_fact* &
542                       1.e12*drygrid(ix,jy)/area(ix,jy)
543!                sparse_dump_u(sp_count_r)=
544!    +            ,1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
545                else ! concentration is zero
[e200b7a]546                  sp_zer=.true.
[b069789]547                endif
548              end do
[e200b7a]549            end do
[b069789]550          else
551            sp_count_i=0
552            sp_count_r=0
553          endif
554          write(unitoutgridppt) sp_count_i
555          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
556          write(unitoutgridppt) sp_count_r
557          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
558!       write(unitoutgridppt) sp_count_u
559!       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
560
561
562! Mixing ratios
563          sp_count_i=0
564          sp_count_r=0
565          sp_fact=-1.
566          sp_zer=.true.
[e200b7a]567          do kz=1,numzgrid
568            do jy=0,numygrid-1
569              do ix=0,numxgrid-1
570                if (grid(ix,jy,kz).gt.smallnum) then
571                  if (sp_zer.eqv..true.) then ! first non zero value
572                    sp_count_i=sp_count_i+1
573                    sparse_dump_i(sp_count_i)= &
574                         ix+jy*numxgrid+kz*numxgrid*numygrid
575                    sp_zer=.false.
576                    sp_fact=sp_fact*(-1.)
[b069789]577                  endif
578                  sp_count_r=sp_count_r+1
579                  sparse_dump_r(sp_count_r)= &
580                       sp_fact* &
581                       1.e12*grid(ix,jy,kz) &
582                       /volume(ix,jy,kz)/outnum* &
583                       weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
584!                sparse_dump_u(sp_count_r)=
585!+              ,1.e12*gridsigma(ix,jy,kz,ks,kp,nage)/volume(ix,jy,kz)/
586!+              outnum*weightair/weightmolar(ks)/
587!+              densityoutgrid(ix,jy,kz)
588                else ! concentration is zero
[e200b7a]589                  sp_zer=.true.
[b069789]590                endif
[e200b7a]591              end do
592            end do
593          end do
[b069789]594          write(unitoutgridppt) sp_count_i
595          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
596          write(unitoutgridppt) sp_count_r
597          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
598!       write(unitoutgridppt) sp_count_u
599!       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
[e200b7a]600
[b069789]601        endif ! output for ppt
[e200b7a]602
[b069789]603      end do
604    end do
[e200b7a]605
606    close(unitoutgridppt)
607    close(unitoutgrid)
608
609  end do
610
[2eefa58]611! RLT Aug 2017
612! Write out conversion factor for dry air
613  inquire(file=path(2)(1:length(2))//'factor_drygrid',exist=lexist)
614  if (lexist) then
615    ! open and append
616    open(unitoutfactor,file=path(2)(1:length(2))//'factor_drygrid',form='unformatted',&
617            status='old',action='write',access='append')
618  else
619    ! create new
620    open(unitoutfactor,file=path(2)(1:length(2))//'factor_drygrid',form='unformatted',&
621            status='new',action='write')
622  endif
623  sp_count_i=0
624  sp_count_r=0
625  sp_fact=-1.
626  sp_zer=.true.
627  do kz=1,numzgrid
628    do jy=0,numygrid-1
629      do ix=0,numxgrid-1
630        if (factor_drygrid(ix,jy,kz).gt.(1.+smallnum).or.factor_drygrid(ix,jy,kz).lt.(1.-smallnum)) then
631          if (sp_zer.eqv..true.) then ! first value not equal to one
632            sp_count_i=sp_count_i+1
633            sparse_dump_i(sp_count_i)= &
634                  ix+jy*numxgrid+kz*numxgrid*numygrid
635            sp_zer=.false.
636            sp_fact=sp_fact*(-1.)
637          endif
638          sp_count_r=sp_count_r+1
639          sparse_dump_r(sp_count_r)= &
640               sp_fact*factor_drygrid(ix,jy,kz)
641        else ! factor is one
642          sp_zer=.true.
643        endif
644      end do
645    end do
646  end do
647  write(unitoutfactor) sp_count_i
648  write(unitoutfactor) (sparse_dump_i(i),i=1,sp_count_i)
649  write(unitoutfactor) sp_count_r
650  write(unitoutfactor) (sparse_dump_r(i),i=1,sp_count_r)
651  close(unitoutfactor)
652
653
[e200b7a]654  if (gridtotal.gt.0.) gridtotalunc=gridsigmatotal/gridtotal
655  if (wetgridtotal.gt.0.) wetgridtotalunc=wetgridsigmatotal/ &
656       wetgridtotal
657  if (drygridtotal.gt.0.) drygridtotalunc=drygridsigmatotal/ &
658       drygridtotal
659
[b069789]660! Dump of receptor concentrations
[e200b7a]661
[b069789]662  if (numreceptor.gt.0 .and. (iout.eq.2 .or. iout.eq.3)  ) then
663    write(unitoutreceptppt) itime
664    do ks=1,nspec
665      write(unitoutreceptppt) (1.e12*creceptor(i,ks)/outnum* &
666           weightair/weightmolar(ks)/densityoutrecept(i),i=1,numreceptor)
667    end do
668  endif
[e200b7a]669
[b069789]670! Dump of receptor concentrations
[e200b7a]671
[b069789]672  if (numreceptor.gt.0) then
673    write(unitoutrecept) itime
674    do ks=1,nspec
675      write(unitoutrecept) (1.e12*creceptor(i,ks)/outnum, &
676           i=1,numreceptor)
677    end do
678  endif
[e200b7a]679
[2eefa58]680! RLT Aug 2017
681! Write out conversion factor for dry air
682  if (numreceptor.gt.0) then
683    inquire(file=path(2)(1:length(2))//'factor_dryreceptor',exist=lexist)
684     if (lexist) then
685     ! open and append
686      open(unitoutfactor,file=path(2)(1:length(2))//'factor_dryreceptor',form='unformatted',&
687              status='old',action='write',access='append')
688    else
689      ! create new
690      open(unitoutfactor,file=path(2)(1:length(2))//'factor_dryreceptor',form='unformatted',&
691              status='new',action='write')
692    endif
693    write(unitoutfactor) itime
694    write(unitoutfactor) (factor_dryrecept(i),i=1,numreceptor)
695    close(unitoutfactor)
696  endif
[e200b7a]697
[b069789]698! Reinitialization of grid
699!*************************
[e200b7a]700
[16b61a5]701! do ks=1,nspec
702!   do kp=1,maxpointspec_act
703!     do i=1,numreceptor
704!       creceptor(i,ks)=0.
705!     end do
706!     do jy=0,numygrid-1
707!       do ix=0,numxgrid-1
708!         do l=1,nclassunc
709!           do nage=1,nageclass
710!             do kz=1,numzgrid
711!               gridunc(ix,jy,kz,ks,kp,l,nage)=0.
712!             end do
713!           end do
714!         end do
715!       end do
716!     end do
717!   end do
718! end do
[db712a8]719  creceptor(:,:)=0.
720  gridunc(:,:,:,:,:,:,:)=0.
[e200b7a]721
722
723end subroutine concoutput
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG