source: flexpart.git/src/concoutput_inversion_nest.f90 @ 2eefa58

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 2eefa58 was 2eefa58, checked in by Espen Sollum ATMOS <eso@…>, 5 years ago

Added Ronas changes for inversion output

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