source: flexpart.git/src/concoutput_surf_nest.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

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