source: flexpart.git/src/concoutput_nest_mpi.f90 @ f3054ea

GFS_025dev
Last change on this file since f3054ea was f3054ea, checked in by Espen Sollum <eso@…>, 4 years ago

Changed from grib_api to eccodes. MPI: implemented linversionout=1; fixed calculation of grid_initial fields.

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