source: flexpart.git/src/concoutput_mpi.f90 @ 8a65cb0

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 8a65cb0 was 8a65cb0, checked in by Espen Sollum ATMOS <espen@…>, 9 years ago

Added code, makefile for dev branch

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