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