source: flexpart.git/src/concoutput_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@…>, 5 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 21.5 KB
RevLine 
[e200b7a]1subroutine concoutput_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
[6a678e3]42  use mean_mod
[e200b7a]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
[2eefa58]51! RLT
52  real :: densitydryrecept(maxreceptor)
53  real :: factor_dryrecept(maxreceptor)
[e200b7a]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)
[6a678e3]73  real(dep_prec) :: auxgrid(nclassunc)
[e200b7a]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
[2eefa58]80  logical :: lexist
[5f9d14a]81  integer :: mind
82! mind        eso:added to ensure identical results between 2&3-fields versions
[e200b7a]83
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
[5f9d14a]122  mind=memind(2)
[e200b7a]123  do kz=1,numzgrid
124    if (kz.eq.1) then
125      halfheight=outheight(1)/2.
126    else
127      halfheight=(outheight(kz)+outheight(kz-1))/2.
128    endif
129    do kzz=2,nz
130      if ((height(kzz-1).lt.halfheight).and. &
131           (height(kzz).gt.halfheight)) goto 46
132    end do
13346   kzz=max(min(kzz,nz),2)
134    dz1=halfheight-height(kzz-1)
135    dz2=height(kzz)-halfheight
136    dz=dz1+dz2
137    do jy=0,numygridn-1
138      do ix=0,numxgridn-1
139        xl=outlon0n+real(ix)*dxoutn
140        yl=outlat0n+real(jy)*dyoutn
141        xl=(xl-xlon0)/dx
142        yl=(yl-ylat0)/dy
143        iix=max(min(nint(xl),nxmin1),0)
144        jjy=max(min(nint(yl),nymin1),0)
[5f9d14a]145        ! densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
146        !      rho(iix,jjy,kzz-1,2)*dz2)/dz
147        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,mind)*dz1+ &
148             rho(iix,jjy,kzz-1,mind)*dz2)/dz
[2eefa58]149! RLT
150        densitydrygrid(ix,jy,kz)=(rho_dry(iix,jjy,kzz,mind)*dz1+ &
151             rho_dry(iix,jjy,kzz-1,mind)*dz2)/dz
[e200b7a]152      end do
153    end do
154  end do
155
[5f9d14a]156  do i=1,numreceptor
157    xl=xreceptor(i)
158    yl=yreceptor(i)
159    iix=max(min(nint(xl),nxmin1),0)
160    jjy=max(min(nint(yl),nymin1),0)
161    !densityoutrecept(i)=rho(iix,jjy,1,2)
162    densityoutrecept(i)=rho(iix,jjy,1,mind)
[2eefa58]163! RLT
164    densitydryrecept(i)=rho_dry(iix,jjy,1,mind)
[5f9d14a]165  end do
[e200b7a]166
[2eefa58]167! RLT
168! conversion factor for output relative to dry air
169  factor_drygrid=densityoutgrid/densitydrygrid
170  factor_dryrecept=densityoutrecept/densitydryrecept
[e200b7a]171
172  ! Output is different for forward and backward simulations
173    do kz=1,numzgrid
174      do jy=0,numygridn-1
175        do ix=0,numxgridn-1
176          if (ldirect.eq.1) then
177            factor3d(ix,jy,kz)=1.e12/volumen(ix,jy,kz)/outnum
178          else
179            factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
180          endif
181        end do
182      end do
183    end do
184
185  !*********************************************************************
186  ! Determine the standard deviation of the mean concentration or mixing
187  ! ratio (uncertainty of the output) and the dry and wet deposition
188  !*********************************************************************
189
190  do ks=1,nspec
191
192  write(anspec,'(i3.3)') ks
[7d02c2f]193 
194  if (DRYBKDEP.or.WETBKDEP) then !scavdep output
195      if (DRYBKDEP) &
196      open(unitoutgrid,file=path(2)(1:length(2))//'grid_drydep_nest_'//adate// &
197           atime//'_'//anspec,form='unformatted')
198      if (WETBKDEP) &
199      open(unitoutgrid,file=path(2)(1:length(2))//'grid_wetdep_nest_'//adate// &
200           atime//'_'//anspec,form='unformatted')
201      write(unitoutgrid) itime
202  else
203     if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
204       if (ldirect.eq.1) then
205         open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_nest_' &
[e200b7a]206           //adate// &
207           atime//'_'//anspec,form='unformatted')
[7d02c2f]208       else
209         open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_nest_' &
[e200b7a]210           //adate// &
211           atime//'_'//anspec,form='unformatted')
[7d02c2f]212       endif
213       write(unitoutgrid) itime
214     endif
215  endif
[e200b7a]216
217  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
218   open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_nest_' &
219        //adate// &
220        atime//'_'//anspec,form='unformatted')
221
222    write(unitoutgridppt) itime
223  endif
224
225  do kp=1,maxpointspec_act
226  do nage=1,nageclass
227
228    do jy=0,numygridn-1
229      do ix=0,numxgridn-1
230
231  ! WET DEPOSITION
232        if ((WETDEP).and.(ldirect.gt.0)) then
233            do l=1,nclassunc
234              auxgrid(l)=wetgriduncn(ix,jy,ks,kp,l,nage)
235            end do
236            call mean(auxgrid,wetgrid(ix,jy), &
237                 wetgridsigma(ix,jy),nclassunc)
238  ! Multiply by number of classes to get total concentration
239            wetgrid(ix,jy)=wetgrid(ix,jy) &
240                 *nclassunc
241  ! Calculate standard deviation of the mean
242            wetgridsigma(ix,jy)= &
243                 wetgridsigma(ix,jy)* &
244                 sqrt(real(nclassunc))
245        endif
246
247  ! DRY DEPOSITION
248        if ((DRYDEP).and.(ldirect.gt.0)) then
249            do l=1,nclassunc
250              auxgrid(l)=drygriduncn(ix,jy,ks,kp,l,nage)
251            end do
252            call mean(auxgrid,drygrid(ix,jy), &
253                 drygridsigma(ix,jy),nclassunc)
254  ! Multiply by number of classes to get total concentration
255            drygrid(ix,jy)=drygrid(ix,jy)* &
256                 nclassunc
257  ! Calculate standard deviation of the mean
258            drygridsigma(ix,jy)= &
259                 drygridsigma(ix,jy)* &
260                 sqrt(real(nclassunc))
261        endif
262
263  ! CONCENTRATION OR MIXING RATIO
264        do kz=1,numzgrid
265            do l=1,nclassunc
266              auxgrid(l)=griduncn(ix,jy,kz,ks,kp,l,nage)
267            end do
268            call mean(auxgrid,grid(ix,jy,kz), &
269                 gridsigma(ix,jy,kz),nclassunc)
270  ! Multiply by number of classes to get total concentration
271            grid(ix,jy,kz)= &
272                 grid(ix,jy,kz)*nclassunc
273  ! Calculate standard deviation of the mean
274            gridsigma(ix,jy,kz)= &
275                 gridsigma(ix,jy,kz)* &
276                 sqrt(real(nclassunc))
277        end do
278      end do
279    end do
280
281
282  !*******************************************************************
283  ! Generate output: may be in concentration (ng/m3) or in mixing
284  ! ratio (ppt) or both
285  ! Output the position and the values alternated multiplied by
286  ! 1 or -1, first line is number of values, number of positions
287  ! For backward simulations, the unit is seconds, stored in grid_time
288  !*******************************************************************
289
290  ! Concentration output
291  !*********************
292  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
293
294  ! Wet deposition
295         sp_count_i=0
296         sp_count_r=0
297         sp_fact=-1.
298         sp_zer=.true.
299         if ((ldirect.eq.1).and.(WETDEP)) then
300         do jy=0,numygridn-1
301            do ix=0,numxgridn-1
302  !oncentraion greater zero
303              if (wetgrid(ix,jy).gt.smallnum) then
304                 if (sp_zer.eqv..true.) then ! first non zero value
305                    sp_count_i=sp_count_i+1
306                    sparse_dump_i(sp_count_i)=ix+jy*numxgridn
307                    sp_zer=.false.
308                    sp_fact=sp_fact*(-1.)
309                 endif
310                 sp_count_r=sp_count_r+1
311                 sparse_dump_r(sp_count_r)= &
312                      sp_fact*1.e12*wetgrid(ix,jy)/arean(ix,jy)
313  !                sparse_dump_u(sp_count_r)=
314  !+                1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
315              else ! concentration is zero
316                  sp_zer=.true.
317              endif
318            end do
319         end do
320         else
321            sp_count_i=0
322            sp_count_r=0
323         endif
324         write(unitoutgrid) sp_count_i
325         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
326         write(unitoutgrid) sp_count_r
327         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
328  !       write(unitoutgrid) sp_count_u
329  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
330
331  ! Dry deposition
332         sp_count_i=0
333         sp_count_r=0
334         sp_fact=-1.
335         sp_zer=.true.
336         if ((ldirect.eq.1).and.(DRYDEP)) then
337          do jy=0,numygridn-1
338            do ix=0,numxgridn-1
339              if (drygrid(ix,jy).gt.smallnum) then
340                 if (sp_zer.eqv..true.) then ! first non zero value
341                    sp_count_i=sp_count_i+1
342                    sparse_dump_i(sp_count_i)=ix+jy*numxgridn
343                    sp_zer=.false.
344                    sp_fact=sp_fact*(-1.)
345                 endif
346                 sp_count_r=sp_count_r+1
347                 sparse_dump_r(sp_count_r)= &
348                      sp_fact* &
349                      1.e12*drygrid(ix,jy)/arean(ix,jy)
350  !                sparse_dump_u(sp_count_r)=
351  !+                1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
352              else ! concentration is zero
353                  sp_zer=.true.
354              endif
355            end do
356          end do
357         else
358            sp_count_i=0
359            sp_count_r=0
360         endif
361         write(unitoutgrid) sp_count_i
362         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
363         write(unitoutgrid) sp_count_r
364         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
365  !       write(*,*) sp_count_u
366  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
367
368
369
370  ! Concentrations
371         sp_count_i=0
372         sp_count_r=0
373         sp_fact=-1.
374         sp_zer=.true.
375          do kz=1,numzgrid
376            do jy=0,numygridn-1
377              do ix=0,numxgridn-1
378                if (grid(ix,jy,kz).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)= &
382                         ix+jy*numxgridn+kz*numxgridn*numygridn
383                    sp_zer=.false.
384                    sp_fact=sp_fact*(-1.)
385                   endif
386                   sp_count_r=sp_count_r+1
387                   sparse_dump_r(sp_count_r)= &
388                        sp_fact* &
389                        grid(ix,jy,kz)* &
390                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
391  !                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
392  !    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
393  !                sparse_dump_u(sp_count_r)=
394  !+               ,gridsigma(ix,jy,kz,ks,kp,nage)*
395  !+               factor(ix,jy,kz)/tot_mu(ks,kp)
396              else ! concentration is zero
397                  sp_zer=.true.
398              endif
399              end do
400            end do
401          end do
402         write(unitoutgrid) sp_count_i
403         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
404         write(unitoutgrid) sp_count_r
405         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
406  !       write(unitoutgrid) sp_count_u
407  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
408
409
410
411    endif !  concentration output
412
413  ! Mixing ratio output
414  !********************
415
416  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
417
418  ! Wet deposition
419         sp_count_i=0
420         sp_count_r=0
421         sp_fact=-1.
422         sp_zer=.true.
423         if ((ldirect.eq.1).and.(WETDEP)) then
424          do jy=0,numygridn-1
425            do ix=0,numxgridn-1
426                if (wetgrid(ix,jy).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
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                      1.e12*wetgrid(ix,jy)/arean(ix,jy)
438  !                sparse_dump_u(sp_count_r)=
439  !    +            ,1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
440              else ! concentration is zero
441                  sp_zer=.true.
442              endif
443            end do
444          end do
445         else
446           sp_count_i=0
447           sp_count_r=0
448         endif
449         write(unitoutgridppt) sp_count_i
450         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
451         write(unitoutgridppt) sp_count_r
452         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
453  !       write(unitoutgridppt) sp_count_u
454  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
455
456
457  ! Dry deposition
458         sp_count_i=0
459         sp_count_r=0
460         sp_fact=-1.
461         sp_zer=.true.
462         if ((ldirect.eq.1).and.(DRYDEP)) then
463          do jy=0,numygridn-1
464            do ix=0,numxgridn-1
465                if (drygrid(ix,jy).gt.smallnum) then
466                  if (sp_zer.eqv..true.) then ! first non zero value
467                    sp_count_i=sp_count_i+1
468                    sparse_dump_i(sp_count_i)= &
469                         ix+jy*numxgridn
470                    sp_zer=.false.
471                    sp_fact=sp_fact*(-1)
472                 endif
473                 sp_count_r=sp_count_r+1
474                 sparse_dump_r(sp_count_r)= &
475                      sp_fact* &
476                      1.e12*drygrid(ix,jy)/arean(ix,jy)
477  !                sparse_dump_u(sp_count_r)=
478  !    +            ,1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
479              else ! concentration is zero
480                  sp_zer=.true.
481              endif
482            end do
483          end do
484         else
485           sp_count_i=0
486           sp_count_r=0
487         endif
488         write(unitoutgridppt) sp_count_i
489         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
490         write(unitoutgridppt) sp_count_r
491         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
492  !       write(unitoutgridppt) sp_count_u
493  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
494
495
496  ! Mixing ratios
497         sp_count_i=0
498         sp_count_r=0
499         sp_fact=-1.
500         sp_zer=.true.
501          do kz=1,numzgrid
502            do jy=0,numygridn-1
503              do ix=0,numxgridn-1
504                if (grid(ix,jy,kz).gt.smallnum) then
505                  if (sp_zer.eqv..true.) then ! first non zero value
506                    sp_count_i=sp_count_i+1
507                    sparse_dump_i(sp_count_i)= &
508                         ix+jy*numxgridn+kz*numxgridn*numygridn
509                    sp_zer=.false.
510                    sp_fact=sp_fact*(-1.)
511                 endif
512                 sp_count_r=sp_count_r+1
513                 sparse_dump_r(sp_count_r)= &
514                      sp_fact* &
515                      1.e12*grid(ix,jy,kz) &
516                      /volumen(ix,jy,kz)/outnum* &
517                      weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
518  !                sparse_dump_u(sp_count_r)=
519  !+              ,1.e12*gridsigma(ix,jy,kz,ks,kp,nage)/volume(ix,jy,kz)/
520  !+              outnum*weightair/weightmolar(ks)/
521  !+              densityoutgrid(ix,jy,kz)
522              else ! concentration is zero
523                  sp_zer=.true.
524              endif
525              end do
526            end do
527          end do
528         write(unitoutgridppt) sp_count_i
529         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
530         write(unitoutgridppt) sp_count_r
531         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
532  !       write(unitoutgridppt) sp_count_u
533  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
534
535      endif ! output for ppt
536
537  end do
538  end do
539
540    close(unitoutgridppt)
541    close(unitoutgrid)
542
543  end do
544
[2eefa58]545! RLT Aug 2017
546! Write out conversion factor for dry air
547  inquire(file=path(2)(1:length(2))//'factor_drygrid_nest',exist=lexist)
548  if (lexist) then
549    ! open and append
550    open(unitoutfactor,file=path(2)(1:length(2))//'factor_drygrid_nest',form='unformatted',&
551            status='old',action='write',access='append')
552  else
553    ! create new
554    open(unitoutfactor,file=path(2)(1:length(2))//'factor_drygrid_nest',form='unformatted',&
555            status='new',action='write')
556  endif
557  sp_count_i=0
558  sp_count_r=0
559  sp_fact=-1.
560  sp_zer=.true.
561  do kz=1,numzgrid
562    do jy=0,numygridn-1
563      do ix=0,numxgridn-1
564        if (factor_drygrid(ix,jy,kz).gt.(1.+smallnum).or.factor_drygrid(ix,jy,kz).lt.(1.-smallnum)) then
565          if (sp_zer.eqv..true.) then ! first value not equal to one
566            sp_count_i=sp_count_i+1
567            sparse_dump_i(sp_count_i)= &
568                  ix+jy*numxgridn+kz*numxgridn*numygridn
569            sp_zer=.false.
570            sp_fact=sp_fact*(-1.)
571          endif
572          sp_count_r=sp_count_r+1
573          sparse_dump_r(sp_count_r)= &
574               sp_fact*factor_drygrid(ix,jy,kz)
575        else ! factor is one
576          sp_zer=.true.
577        endif
578      end do
579    end do
580  end do
581  write(unitoutfactor) sp_count_i
582  write(unitoutfactor) (sparse_dump_i(i),i=1,sp_count_i)
583  write(unitoutfactor) sp_count_r
584  write(unitoutfactor) (sparse_dump_r(i),i=1,sp_count_r)
585  close(unitoutfactor)
[e200b7a]586
587
588  ! Reinitialization of grid
589  !*************************
590
[db712a8]591  ! do ks=1,nspec
592  ! do kp=1,maxpointspec_act
593  !   do i=1,numreceptor
594  !     creceptor(i,ks)=0.
595  !   end do
596  !   do jy=0,numygridn-1
597  !     do ix=0,numxgridn-1
598  !       do l=1,nclassunc
599  !         do nage=1,nageclass
600  !           do kz=1,numzgrid
601  !             griduncn(ix,jy,kz,ks,kp,l,nage)=0.
602  !           end do
603  !         end do
604  !       end do
605  !     end do
606  !   end do
607  ! end do
608  ! end do
609  creceptor(:,:)=0.
610  griduncn(:,:,:,:,:,:,:)=0.
[e200b7a]611
612
613end subroutine concoutput_nest
614
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG