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