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