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

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

Merge branch 'dev' into espen

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