source: flexpart.git/src/concoutput_inversion.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: 29.6 KB
Line 
1subroutine concoutput_inversion(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!     January 2017,  Separate files by release but include all timesteps
30!                                                                            *
31!*****************************************************************************
32!                                                                            *
33! Variables:                                                                 *
34! outnum          number of samples                                          *
35! ncells          number of cells with non-zero concentrations               *
36! sparse          .true. if in sparse matrix format, else .false.            *
37! tot_mu          1 for forward, initial mass mixing ration for backw. runs  *
38!                                                                            *
39!*****************************************************************************
40
41  use unc_mod
42  use point_mod
43  use outg_mod
44  use par_mod
45  use com_mod
46  use mean_mod
47
48  implicit none
49
50  real(kind=dp) :: jul
51  integer :: itime,i,ix,jy,kz,ks,kp,l,iix,jjy,kzz,nage,jjjjmmdd,ihmmss
52  integer :: sp_count_i,sp_count_r
53  real :: sp_fact
54  real :: outnum,densityoutrecept(maxreceptor),xl,yl
55! RLT
56  real :: densitydryrecept(maxreceptor)
57  real :: factor_dryrecept(maxreceptor)
58
59!real densityoutgrid(0:numxgrid-1,0:numygrid-1,numzgrid),
60!    +grid(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,maxpointspec_act,
61!    +    maxageclass)
62!real wetgrid(0:numxgrid-1,0:numygrid-1,maxspec,maxpointspec_act,
63!    +       maxageclass)
64!real drygrid(0:numxgrid-1,0:numygrid-1,maxspec,
65!    +       maxpointspec_act,maxageclass)
66!real gridsigma(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,
67!    +       maxpointspec_act,maxageclass),
68!    +     drygridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
69!    +     maxpointspec_act,maxageclass),
70!    +     wetgridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
71!    +     maxpointspec_act,maxageclass)
72!real factor(0:numxgrid-1,0:numygrid-1,numzgrid)
73!real sparse_dump_r(numxgrid*numygrid*numzgrid)
74!integer sparse_dump_i(numxgrid*numygrid*numzgrid)
75!real sparse_dump_u(numxgrid*numygrid*numzgrid)
76
77  real(dep_prec) :: auxgrid(nclassunc)
78  real(sp) :: gridtotal,gridsigmatotal,gridtotalunc
79  real(dep_prec) :: wetgridtotal,wetgridsigmatotal,wetgridtotalunc
80  real(dep_prec) :: drygridtotal,drygridsigmatotal,drygridtotalunc
81  real :: halfheight,dz,dz1,dz2,tot_mu(maxspec,maxpointspec_act)
82  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
83  real,parameter :: weightair=28.97
84  logical :: sp_zer
85  character :: adate*8,atime*6
86  character(len=3) :: anspec
87  logical :: lexist
88  character :: areldate*8,areltime*6
89  logical,save :: lstart=.true.
90  logical,save,allocatable,dimension(:) :: lstartrel
91  integer :: ierr
92  character(LEN=100) :: dates_char
93  integer, parameter :: unitrelnames=654
94
95
96  if(lstart) then
97    allocate(lstartrel(maxpointspec_act))
98    lstartrel(:)=.true.
99  endif
100  print*, 'lstartrel = ',lstartrel
101
102  if (verbosity.eq.1) then
103    print*,'inside concoutput_inversion '
104    CALL SYSTEM_CLOCK(count_clock)
105    WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
106  endif
107
108! Determine current calendar date
109!**********************************************************
110
111  jul=bdate+real(itime,kind=dp)/86400._dp
112  call caldate(jul,jjjjmmdd,ihmmss)
113  write(adate,'(i8.8)') jjjjmmdd
114  write(atime,'(i6.6)') ihmmss
115!  write(unitdates,'(a)') adate//atime
116
117
118! Prepare output files for dates
119!**********************************************************
120
121  ! Overwrite existing dates file on first call, later append to it
122  ! If 'dates' file exists in output directory, copy to new file dates.old
123  inquire(file=path(2)(1:length(2))//'dates', exist=lexist)
124  if (lexist.and.lstart) then
125    ! copy contents of existing dates file to dates.old
126    print*, 'warning: replacing old dates file'
127    open(unit=unitdates, file=path(2)(1:length(2))//'dates',form='formatted', &
128         &access='sequential', status='old', action='read', iostat=ierr)
129    open(unit=unittmp, file=path(2)(1:length(2))//'dates.old', access='sequential', &
130         &status='replace', action='write', form='formatted', iostat=ierr)
131    do while (.true.)
132      read(unitdates, '(a)', iostat=ierr) dates_char
133      if (ierr.ne.0) exit
134      write(unit=unittmp, fmt='(a)', iostat=ierr, advance='yes') trim(dates_char)
135    end do
136    close(unit=unitdates)
137    close(unit=unittmp)
138    ! create new dates file
139    open(unit=unitdates, file=path(2)(1:length(2))//'dates',form='formatted', &
140         &access='sequential', status='replace', iostat=ierr)
141    close(unit=unitdates)
142  endif
143
144  open(unitdates,file=path(2)(1:length(2))//'dates', ACCESS='APPEND')
145  write(unitdates,'(a)') adate//atime
146  close(unitdates)
147
148  !CGZ: Make a filename with names of releases
149  if (lstart) then
150    open(unit=unitrelnames, file=path(2)(1:length(2))//'releases_out',form='formatted', &
151         &access='sequential', status='replace', iostat=ierr)
152    close(unitrelnames)
153  endif
154
155  print*, 'after creating dates files: lstart = ',lstart
156!  print*, 'outnum:',outnum
157!  print*, 'datetime:',adate//atime
158
159
160! For forward simulations, output fields have dimension MAXSPEC,
161! for backward simulations, output fields have dimension MAXPOINT.
162! Thus, make loops either about nspec, or about numpoint
163!*****************************************************************
164
165
166  if (ldirect.eq.1) then
167    do ks=1,nspec
168      do kp=1,maxpointspec_act
169        tot_mu(ks,kp)=1
170      end do
171    end do
172  else
173    do ks=1,nspec
174      do kp=1,maxpointspec_act
175        tot_mu(ks,kp)=xmass(kp,ks)
176      end do
177    end do
178  endif
179
180
181  if (verbosity.eq.1) then
182    print*,'concoutput_inversion 2'
183    CALL SYSTEM_CLOCK(count_clock)
184    WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
185  endif
186
187!*******************************************************************
188! Compute air density: sufficiently accurate to take it
189! from coarse grid at some time
190! Determine center altitude of output layer, and interpolate density
191! data to that altitude
192!*******************************************************************
193
194  do kz=1,numzgrid
195    if (kz.eq.1) then
196      halfheight=outheight(1)/2.
197    else
198      halfheight=(outheight(kz)+outheight(kz-1))/2.
199    endif
200    do kzz=2,nz
201      if ((height(kzz-1).lt.halfheight).and. &
202           (height(kzz).gt.halfheight)) goto 46
203    end do
20446  kzz=max(min(kzz,nz),2)
205    dz1=halfheight-height(kzz-1)
206    dz2=height(kzz)-halfheight
207    dz=dz1+dz2
208    do jy=0,numygrid-1
209      do ix=0,numxgrid-1
210        xl=outlon0+real(ix)*dxout
211        yl=outlat0+real(jy)*dyout
212        xl=(xl-xlon0)/dx
213        yl=(yl-ylat0)/dy
214        iix=max(min(nint(xl),nxmin1),0)
215        jjy=max(min(nint(yl),nymin1),0)
216        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
217             rho(iix,jjy,kzz-1,2)*dz2)/dz
218! RLT
219        densitydrygrid(ix,jy,kz)=(rho_dry(iix,jjy,kzz,2)*dz1+ &
220             rho_dry(iix,jjy,kzz-1,2)*dz2)/dz
221      end do
222    end do
223  end do
224
225  do i=1,numreceptor
226    xl=xreceptor(i)
227    yl=yreceptor(i)
228    iix=max(min(nint(xl),nxmin1),0)
229    jjy=max(min(nint(yl),nymin1),0)
230    densityoutrecept(i)=rho(iix,jjy,1,2)
231! RLT
232    densitydryrecept(i)=rho_dry(iix,jjy,1,2)
233  end do
234
235! RLT
236! conversion factor for output relative to dry air
237  factor_drygrid=densityoutgrid/densitydrygrid
238  factor_dryrecept=densityoutrecept/densitydryrecept
239
240! Output is different for forward and backward simulations
241  do kz=1,numzgrid
242    do jy=0,numygrid-1
243      do ix=0,numxgrid-1
244        if (ldirect.eq.1) then
245          factor3d(ix,jy,kz)=1.e12/volume(ix,jy,kz)/outnum
246        else
247          factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
248        endif
249      end do
250    end do
251  end do
252
253!*********************************************************************
254! Determine the standard deviation of the mean concentration or mixing
255! ratio (uncertainty of the output) and the dry and wet deposition
256!*********************************************************************
257
258  if (verbosity.eq.1) then
259    print*,'concoutput_inversion 3 (sd)'
260    CALL SYSTEM_CLOCK(count_clock)
261    WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
262  endif
263  gridtotal=0.
264  gridsigmatotal=0.
265  gridtotalunc=0.
266  wetgridtotal=0.
267  wetgridsigmatotal=0.
268  wetgridtotalunc=0.
269  drygridtotal=0.
270  drygridsigmatotal=0.
271  drygridtotalunc=0.
272
273  do ks=1,nspec
274
275    write(anspec,'(i3.3)') ks
276
277    ! loop over releases
278    do kp=1,maxpointspec_act
279
280      print*, 'itime = ',itime
281      !print*, 'lage(1) = ',lage(1)
282      print*, 'ireleasestart(kp) = ',ireleasestart(kp)
283      print*, 'ireleaseend(kp) = ',ireleaseend(kp)
284
285      ! check itime is within release and backward trajectory length
286      if (nageclass.eq.1) then
287        if ((itime.gt.ireleaseend(kp)).or.(itime.lt.(ireleasestart(kp)-lage(1)))) then
288          go to 10
289        endif
290      endif
291
292      ! calculate date of release for filename
293      jul=bdate+real(ireleasestart(kp),kind=dp)/86400._dp    ! this is the current day
294      call caldate(jul,jjjjmmdd,ihmmss)     
295      write(areldate,'(i8.8)') jjjjmmdd
296      write(areltime,'(i6.6)') ihmmss
297      print*, 'areldate/areltime = ',areldate//areltime
298
299      ! calculate date of field
300      jul=bdate+real(itime,kind=dp)/86400._dp
301      call caldate(jul,jjjjmmdd,ihmmss)
302      write(adate,'(i8.8)') jjjjmmdd
303      write(atime,'(i6.6)') ihmmss
304!      print*, adate//atime
305
306      if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
307        if (ldirect.eq.1) then
308          ! concentrations
309          inquire(file=path(2)(1:length(2))//'grid_conc_'//areldate// &
310                  areltime//'_'//anspec,exist=lexist)
311          if(lexist.and..not.lstartrel(kp)) then
312            ! open and append to existing file
313            open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_'//areldate// &
314                 areltime//'_'//anspec,form='unformatted',status='old',action='write',access='append')
315          else
316            ! open new file
317            open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_'//areldate// &
318                 areltime//'_'//anspec,form='unformatted',status='replace',action='write')
319          endif
320        else
321          ! residence times
322          inquire(file=path(2)(1:length(2))//'grid_time_'//areldate// &
323                  areltime//'_'//anspec,exist=lexist)
324          if(lexist.and..not.lstartrel(kp)) then
325            ! open and append to existing file
326            open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_'//areldate// &
327                 areltime//'_'//anspec,form='unformatted',status='old',action='write',access='append')
328          else
329            ! open new file
330            open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_'//areldate// &
331                 areltime//'_'//anspec,form='unformatted',status='replace',action='write')
332            ! add part of the filename to a file (similar to dates) for easier post-processing
333            open(unit=unitrelnames, file=path(2)(1:length(2))//'releases_out',form='formatted', &
334                 & access='APPEND', iostat=ierr)
335            write(unitrelnames,'(a)') areldate//areltime//'_'//anspec
336            close(unitrelnames)
337          endif
338        endif
339        write(unitoutgrid) jjjjmmdd
340        write(unitoutgrid) ihmmss
341      endif
342
343      if ((iout.eq.2).or.(iout.eq.3)) then     
344        ! mixing ratio
345        inquire(file=path(2)(1:length(2))//'grid_pptv_'//areldate// &
346                areltime//'_'//anspec,exist=lexist)
347        if(lexist.and..not.lstartrel(kp)) then
348          ! open and append to existing file
349          open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_'//areldate// &
350               areltime//'_'//anspec,form='unformatted',status='old',action='write',access='append')
351        else
352          ! open new file
353          open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_'//areldate// &
354               areltime//'_'//anspec,form='unformatted',status='replace',action='write')
355        endif   
356        write(unitoutgridppt) jjjjmmdd
357        write(unitoutgridppt) ihmmss
358      endif
359
360      lstartrel(kp)=.false.
361
362      do nage=1,nageclass
363
364        do jy=0,numygrid-1
365          do ix=0,numxgrid-1
366
367!! WET DEPOSITION
368!            if ((WETDEP).and.(ldirect.gt.0)) then
369!              do l=1,nclassunc
370!                auxgrid(l)=wetgridunc(ix,jy,ks,kp,l,nage)
371!              end do
372!              call mean(auxgrid,wetgrid(ix,jy), &
373!                   wetgridsigma(ix,jy),nclassunc)
374!! Multiply by number of classes to get total concentration
375!              wetgrid(ix,jy)=wetgrid(ix,jy) &
376!                   *nclassunc
377!              wetgridtotal=wetgridtotal+wetgrid(ix,jy)
378!! Calculate standard deviation of the mean
379!              wetgridsigma(ix,jy)= &
380!                   wetgridsigma(ix,jy)* &
381!                   sqrt(real(nclassunc))
382!              wetgridsigmatotal=wetgridsigmatotal+ &
383!                   wetgridsigma(ix,jy)
384!            endif
385
386!! DRY DEPOSITION
387!            if ((DRYDEP).and.(ldirect.gt.0)) then
388!              do l=1,nclassunc
389!                auxgrid(l)=drygridunc(ix,jy,ks,kp,l,nage)
390!              end do
391!              call mean(auxgrid,drygrid(ix,jy), &
392!                   drygridsigma(ix,jy),nclassunc)
393!! Multiply by number of classes to get total concentration
394!              drygrid(ix,jy)=drygrid(ix,jy)* &
395!                   nclassunc
396!              drygridtotal=drygridtotal+drygrid(ix,jy)
397!! Calculate standard deviation of the mean
398!              drygridsigma(ix,jy)= &
399!                   drygridsigma(ix,jy)* &
400!                   sqrt(real(nclassunc))
401!125           drygridsigmatotal=drygridsigmatotal+ &
402!                   drygridsigma(ix,jy)
403!            endif
404
405! CONCENTRATION OR MIXING RATIO
406            do kz=1,numzgrid
407              do l=1,nclassunc
408                auxgrid(l)=gridunc(ix,jy,kz,ks,kp,l,nage)
409              end do
410              call mean(auxgrid,grid(ix,jy,kz), &
411                   gridsigma(ix,jy,kz),nclassunc)
412! Multiply by number of classes to get total concentration
413              grid(ix,jy,kz)= &
414                   grid(ix,jy,kz)*nclassunc
415              gridtotal=gridtotal+grid(ix,jy,kz)
416! Calculate standard deviation of the mean
417              gridsigma(ix,jy,kz)= &
418                   gridsigma(ix,jy,kz)* &
419                   sqrt(real(nclassunc))
420              gridsigmatotal=gridsigmatotal+ &
421                   gridsigma(ix,jy,kz)
422            end do
423          end do
424        end do
425
426
427!*******************************************************************
428! Generate output: may be in concentration (ng/m3) or in mixing
429! ratio (ppt) or both
430! Output the position and the values alternated multiplied by
431! 1 or -1, first line is number of values, number of positions
432! For backward simulations, the unit is seconds, stored in grid_time
433!*******************************************************************
434
435        if (verbosity.eq.1) then
436          print*,'concoutput_inversion 4 (output)'
437          CALL SYSTEM_CLOCK(count_clock)
438          WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
439        endif
440
441! Concentration output
442!*********************
443
444        if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
445
446          if (verbosity.eq.1) then
447            print*,'concoutput_inversion (Wet deposition)'
448            CALL SYSTEM_CLOCK(count_clock)
449            WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
450          endif
451
452! Wet deposition
453!          sp_count_i=0
454!          sp_count_r=0
455!          sp_fact=-1.
456!          sp_zer=.true.
457!          if ((ldirect.eq.1).and.(WETDEP)) then
458!            do jy=0,numygrid-1
459!              do ix=0,numxgrid-1
460!! concentration greater zero
461!                if (wetgrid(ix,jy).gt.smallnum) then
462!                  if (sp_zer.eqv..true.) then ! first non zero value
463!                    sp_count_i=sp_count_i+1
464!                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
465!                    sp_zer=.false.
466!                    sp_fact=sp_fact*(-1.)
467!                  endif
468!                  sp_count_r=sp_count_r+1
469!                  sparse_dump_r(sp_count_r)= &
470!                       sp_fact*1.e12*wetgrid(ix,jy)/area(ix,jy)
471!                  sparse_dump_u(sp_count_r)= &
472!                       1.e12*wetgridsigma(ix,jy)/area(ix,jy)
473!                else ! concentration is zero
474!                  sp_zer=.true.
475!                endif
476!              end do
477!            end do
478!          else
479!            sp_count_i=0
480!            sp_count_r=0
481!          endif
482!          write(unitoutgrid) sp_count_i
483!          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
484!          write(unitoutgrid) sp_count_r
485!          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
486!!         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
487!
488!          if (verbosity.eq.1) then
489!            print*,'concoutput_surf (Dry deposition)'
490!            CALL SYSTEM_CLOCK(count_clock)
491!            WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
492!          endif
493!! Dry deposition
494!          sp_count_i=0
495!          sp_count_r=0
496!          sp_fact=-1.
497!          sp_zer=.true.
498!          if ((ldirect.eq.1).and.(DRYDEP)) then
499!            do jy=0,numygrid-1
500!              do ix=0,numxgrid-1
501!                if (drygrid(ix,jy).gt.smallnum) then
502!                  if (sp_zer.eqv..true.) then ! first non zero value
503!                    sp_count_i=sp_count_i+1
504!                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
505!                    sp_zer=.false.
506!                    sp_fact=sp_fact*(-1.)
507!                  endif
508!                  sp_count_r=sp_count_r+1
509!                  sparse_dump_r(sp_count_r)= &
510!                       sp_fact* &
511!                       1.e12*drygrid(ix,jy)/area(ix,jy)
512!                  sparse_dump_u(sp_count_r)= &
513!                       1.e12*drygridsigma(ix,jy)/area(ix,jy)
514!                else ! concentration is zero
515!                  sp_zer=.true.
516!                endif
517!              end do
518!            end do
519!          else
520!            sp_count_i=0
521!            sp_count_r=0
522!          endif
523!          write(unitoutgrid) sp_count_i
524!          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
525!          write(unitoutgrid) sp_count_r
526!          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
527!!         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
528
529          if (verbosity.eq.1) then
530            print*,'concoutput_inversion (Concentrations)'
531            CALL SYSTEM_CLOCK(count_clock)
532            WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
533          endif
534
535! Concentrations
536
537! surf_only write only 1st layer
538
539          sp_count_i=0
540          sp_count_r=0
541          sp_fact=-1.
542          sp_zer=.true.
543          do kz=1,1
544            do jy=0,numygrid-1
545              do ix=0,numxgrid-1
546                if (grid(ix,jy,kz).gt.smallnum) then
547                  if (sp_zer.eqv..true.) then ! first non zero value
548                    sp_count_i=sp_count_i+1
549                    sparse_dump_i(sp_count_i)= &
550                         ix+jy*numxgrid+kz*numxgrid*numygrid
551                    sp_zer=.false.
552                    sp_fact=sp_fact*(-1.)
553                  endif
554                  sp_count_r=sp_count_r+1
555                  sparse_dump_r(sp_count_r)= &
556                       sp_fact* &
557                       grid(ix,jy,kz)* &
558                       factor3d(ix,jy,kz)/tot_mu(ks,kp)
559!                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
560!    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
561                  sparse_dump_u(sp_count_r)= &
562                       gridsigma(ix,jy,kz)* &
563                       factor3d(ix,jy,kz)/tot_mu(ks,kp)
564
565                else ! concentration is zero
566                  sp_zer=.true.
567                endif
568              end do
569            end do
570          end do
571          write(unitoutgrid) sp_count_i
572          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
573          write(unitoutgrid) sp_count_r
574          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
575!          write(unitoutgrid) sp_count_r
576!         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
577
578        endif !  concentration output
579
580! Mixing ratio output
581!********************
582
583        if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
584!
585!! Wet deposition
586!          sp_count_i=0
587!          sp_count_r=0
588!          sp_fact=-1.
589!          sp_zer=.true.
590!          if ((ldirect.eq.1).and.(WETDEP)) then
591!            do jy=0,numygrid-1
592!              do ix=0,numxgrid-1
593!                if (wetgrid(ix,jy).gt.smallnum) then
594!                  if (sp_zer.eqv..true.) then ! first non zero value
595!                    sp_count_i=sp_count_i+1
596!                    sparse_dump_i(sp_count_i)= &
597!                         ix+jy*numxgrid
598!                    sp_zer=.false.
599!                    sp_fact=sp_fact*(-1.)
600!                  endif
601!                  sp_count_r=sp_count_r+1
602!                  sparse_dump_r(sp_count_r)= &
603!                       sp_fact* &
604!                       1.e12*wetgrid(ix,jy)/area(ix,jy)
605!                  sparse_dump_u(sp_count_r)= &
606!                       1.e12*wetgridsigma(ix,jy)/area(ix,jy)
607!                else ! concentration is zero
608!                  sp_zer=.true.
609!                endif
610!              end do
611!            end do
612!          else
613!            sp_count_i=0
614!            sp_count_r=0
615!          endif
616!          write(unitoutgridppt) sp_count_i
617!          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
618!          write(unitoutgridppt) sp_count_r
619!          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
620!!         write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
621!
622
623! Dry deposition
624!          sp_count_i=0
625!          sp_count_r=0
626!          sp_fact=-1.
627!          sp_zer=.true.
628!          if ((ldirect.eq.1).and.(DRYDEP)) then
629!            do jy=0,numygrid-1
630!              do ix=0,numxgrid-1
631!                if (drygrid(ix,jy).gt.smallnum) then
632!                  if (sp_zer.eqv..true.) then ! first non zero value
633!                    sp_count_i=sp_count_i+1
634!                    sparse_dump_i(sp_count_i)= &
635!                         ix+jy*numxgrid
636!                    sp_zer=.false.
637!                    sp_fact=sp_fact*(-1)
638!                  endif
639!                  sp_count_r=sp_count_r+1
640!                  sparse_dump_r(sp_count_r)= &
641!                       sp_fact* &
642!                       1.e12*drygrid(ix,jy)/area(ix,jy)
643!                  sparse_dump_u(sp_count_r)= &
644!                       1.e12*drygridsigma(ix,jy)/area(ix,jy)
645!                else ! concentration is zero
646!                  sp_zer=.true.
647!                endif
648!              end do
649!            end do
650!          else
651!            sp_count_i=0
652!            sp_count_r=0
653!          endif
654!          write(unitoutgridppt) sp_count_i
655!          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
656!          write(unitoutgridppt) sp_count_r
657!          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
658!!         write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
659!
660
661! Mixing ratios
662
663! surf_only write only 1st layer
664
665          sp_count_i=0
666          sp_count_r=0
667          sp_fact=-1.
668          sp_zer=.true.
669          do kz=1,1
670            do jy=0,numygrid-1
671              do ix=0,numxgrid-1
672                if (grid(ix,jy,kz).gt.smallnum) then
673                  if (sp_zer.eqv..true.) then ! first non zero value
674                    sp_count_i=sp_count_i+1
675                    sparse_dump_i(sp_count_i)= &
676                         ix+jy*numxgrid+kz*numxgrid*numygrid
677                    sp_zer=.false.
678                    sp_fact=sp_fact*(-1.)
679                  endif
680                  sp_count_r=sp_count_r+1
681                  sparse_dump_r(sp_count_r)= &
682                       sp_fact* &
683                       1.e12*grid(ix,jy,kz) &
684                       /volume(ix,jy,kz)/outnum* &
685                       weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
686                  sparse_dump_u(sp_count_r)= &
687                       1.e12*gridsigma(ix,jy,kz)/volume(ix,jy,kz)/ &
688                       outnum*weightair/weightmolar(ks)/ &
689                       densityoutgrid(ix,jy,kz)
690                else ! concentration is zero
691                  sp_zer=.true.
692                endif
693              end do
694            end do
695          end do
696          write(unitoutgridppt) sp_count_i
697          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
698          write(unitoutgridppt) sp_count_r
699          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
700!          write(unitoutgridppt) sp_count_r
701!         write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
702
703        endif ! output for ppt
704
705      end do  ! nageclass
706
707      close(unitoutgridppt)
708      close(unitoutgrid)
709
710      ! itime is outside range
71110    continue
712
713    end do  ! maxpointspec_act
714
715  end do  ! nspec
716
717! RLT Aug 2017
718! Write out conversion factor for dry air
719  inquire(file=path(2)(1:length(2))//'factor_drygrid',exist=lexist)
720  if (lexist.and..not.lstart) then
721    ! open and append
722    open(unitoutfactor,file=path(2)(1:length(2))//'factor_drygrid',form='unformatted',&
723            status='old',action='write',access='append')
724  else
725    ! create new
726    open(unitoutfactor,file=path(2)(1:length(2))//'factor_drygrid',form='unformatted',&
727            status='replace',action='write')
728  endif
729  sp_count_i=0
730  sp_count_r=0
731  sp_fact=-1.
732  sp_zer=.true.
733  do kz=1,1
734    do jy=0,numygrid-1
735      do ix=0,numxgrid-1
736        if (factor_drygrid(ix,jy,kz).gt.(1.+smallnum).or.factor_drygrid(ix,jy,kz).lt.(1.-smallnum)) then
737          if (sp_zer.eqv..true.) then ! first value not equal to one
738            sp_count_i=sp_count_i+1
739            sparse_dump_i(sp_count_i)= &
740                  ix+jy*numxgrid+kz*numxgrid*numygrid
741            sp_zer=.false.
742            sp_fact=sp_fact*(-1.)
743          endif
744          sp_count_r=sp_count_r+1
745          sparse_dump_r(sp_count_r)= &
746               sp_fact*factor_drygrid(ix,jy,kz)
747        else ! factor is one
748          sp_zer=.true.
749        endif
750      end do
751    end do
752  end do
753  write(unitoutfactor) sp_count_i
754  write(unitoutfactor) (sparse_dump_i(i),i=1,sp_count_i)
755  write(unitoutfactor) sp_count_r
756  write(unitoutfactor) (sparse_dump_r(i),i=1,sp_count_r)
757  close(unitoutfactor)
758
759
760  if (gridtotal.gt.0.) gridtotalunc=gridsigmatotal/gridtotal
761!  if (wetgridtotal.gt.0.) wetgridtotalunc=wetgridsigmatotal/ &
762!       wetgridtotal
763!  if (drygridtotal.gt.0.) drygridtotalunc=drygridsigmatotal/ &
764!       drygridtotal
765
766! Dump of receptor concentrations
767
768  if (numreceptor.gt.0 .and. (iout.eq.2 .or. iout.eq.3)  ) then
769    write(unitoutreceptppt) itime
770    do ks=1,nspec
771      write(unitoutreceptppt) (1.e12*creceptor(i,ks)/outnum* &
772           weightair/weightmolar(ks)/densityoutrecept(i),i=1,numreceptor)
773    end do
774  endif
775
776! Dump of receptor concentrations
777
778  if (numreceptor.gt.0) then
779    write(unitoutrecept) itime
780    do ks=1,nspec
781      write(unitoutrecept) (1.e12*creceptor(i,ks)/outnum, &
782           i=1,numreceptor)
783    end do
784  endif
785
786! RLT Aug 2017
787! Write out conversion factor for dry air
788  if (numreceptor.gt.0) then
789    inquire(file=path(2)(1:length(2))//'factor_dryreceptor',exist=lexist)
790     if (lexist.and..not.lstart) then
791     ! open and append
792      open(unitoutfactor,file=path(2)(1:length(2))//'factor_dryreceptor',form='unformatted',&
793              status='old',action='write',access='append')
794    else
795      ! create new
796      open(unitoutfactor,file=path(2)(1:length(2))//'factor_dryreceptor',form='unformatted',&
797              status='replace',action='write')
798    endif
799    write(unitoutfactor) itime
800    write(unitoutfactor) (factor_dryrecept(i),i=1,numreceptor)
801    close(unitoutfactor)
802  endif
803
804  ! reset lstart
805  if (lstart) then
806    lstart=.false.
807  endif
808  print*, 'after writing output files: lstart = ',lstart
809
810
811! Reinitialization of grid
812!*************************
813
814  do ks=1,nspec
815    do kp=1,maxpointspec_act
816      do i=1,numreceptor
817        creceptor(i,ks)=0.
818      end do
819      do jy=0,numygrid-1
820        do ix=0,numxgrid-1
821          do l=1,nclassunc
822            do nage=1,nageclass
823              do kz=1,numzgrid
824                gridunc(ix,jy,kz,ks,kp,l,nage)=0.
825              end do
826            end do
827          end do
828        end do
829      end do
830    end do
831  end do
832
833end subroutine concoutput_inversion
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG