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