source: flexpart.git/src/concoutput.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@…>, 5 years ago

move license from headers to a different file

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