CTBTO FLEXPART WO4 (2015-10-15)
 All Classes Files Functions Variables
concoutput.F90
Go to the documentation of this file.
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 
22 subroutine 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  !*****************************************************************************
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 
66  implicit none
67 
68  real(kind=dp) :: jul
69  integer :: itime,i,ix,jy,kz,ks,kp,l,iix,jjy,kzz,nage,jjjjmmdd,ihmmss
70  integer :: sp_count_i,sp_count_r
71  real :: sp_fact
72  real :: outnum,densityoutrecept(maxreceptor),xl,yl
73 
74  !real densityoutgrid(0:numxgrid-1,0:numygrid-1,numzgrid),
75  ! +grid(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,maxpointspec_act,
76  ! + maxageclass)
77  !real wetgrid(0:numxgrid-1,0:numygrid-1,maxspec,maxpointspec_act,
78  ! + maxageclass)
79  !real drygrid(0:numxgrid-1,0:numygrid-1,maxspec,
80  ! + maxpointspec_act,maxageclass)
81  !real gridsigma(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,
82  ! + maxpointspec_act,maxageclass),
83  ! + drygridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
84  ! + maxpointspec_act,maxageclass),
85  ! + wetgridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
86  ! + maxpointspec_act,maxageclass)
87  !real factor(0:numxgrid-1,0:numygrid-1,numzgrid)
88  !real sparse_dump_r(numxgrid*numygrid*numzgrid)
89  !integer sparse_dump_i(numxgrid*numygrid*numzgrid)
90 
91  !real sparse_dump_u(numxgrid*numygrid*numzgrid)
92  real :: auxgrid(nclassunc),gridtotal,gridsigmatotal,gridtotalunc
93  real :: wetgridtotal,wetgridsigmatotal,wetgridtotalunc
94  real :: drygridtotal,drygridsigmatotal,drygridtotalunc
95  real :: halfheight,dz,dz1,dz2,tot_mu(maxspec,maxpointspec_act)
96  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
97  real,parameter :: weightair=28.97
98  logical :: sp_zer
99  character :: adate*8,atime*6
100  character(len=3) :: anspec
101 
102 
103  ! Determine current calendar date, needed for the file name
104  !**********************************************************
105 
106  jul=bdate+real(itime,kind=dp)/86400._dp
107  call caldate(jul,jjjjmmdd,ihmmss)
108  write(adate,'(i8.8)') jjjjmmdd
109  write(atime,'(i6.6)') ihmmss
110  write(unitdates,'(a)') adate//atime
111 
112 
113  ! For forward simulations, output fields have dimension MAXSPEC,
114  ! for backward simulations, output fields have dimension MAXPOINT.
115  ! Thus, make loops either about nspec, or about numpoint
116  !*****************************************************************
117 
118 
119  if (ldirect.eq.1) then
120  do ks=1,nspec
121  do kp=1,maxpointspec_act
122  tot_mu(ks,kp)=1
123  end do
124  end do
125  else
126  do ks=1,nspec
127  do kp=1,maxpointspec_act
128 #if defined WITH_CTBTO_PATCHES
129  tot_mu(ks,kp)=1
130 #else
131  tot_mu(ks,kp)=xmass(kp,ks)
132 #endif
133  end do
134  end do
135  endif
136 
137 
138  !*******************************************************************
139  ! Compute air density: sufficiently accurate to take it
140  ! from coarse grid at some time
141  ! Determine center altitude of output layer, and interpolate density
142  ! data to that altitude
143  !*******************************************************************
144 
145  do kz=1,numzgrid
146  if (kz.eq.1) then
147  halfheight=outheight(1)/2.
148  else
149  halfheight=(outheight(kz)+outheight(kz-1))/2.
150  endif
151  do kzz=2,nz
152  if ((height(kzz-1).lt.halfheight).and. &
153  (height(kzz).gt.halfheight)) goto 46
154  end do
155 46 kzz=max(min(kzz,nz),2)
156  dz1=halfheight-height(kzz-1)
157  dz2=height(kzz)-halfheight
158  dz=dz1+dz2
159  do jy=0,numygrid-1
160  do ix=0,numxgrid-1
161  xl=outlon0+real(ix)*dxout
162  yl=outlat0+real(jy)*dyout
163  xl=(xl-xlon0)/dx
164  yl=(yl-ylat0)/dx
165  iix=max(min(nint(xl),nxmin1),0)
166  jjy=max(min(nint(yl),nymin1),0)
167  densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
168  rho(iix,jjy,kzz-1,2)*dz2)/dz
169  end do
170  end do
171  end do
172 
173  do i=1,numreceptor
174  xl=xreceptor(i)
175  yl=yreceptor(i)
176  iix=max(min(nint(xl),nxmin1),0)
177  jjy=max(min(nint(yl),nymin1),0)
178  densityoutrecept(i)=rho(iix,jjy,1,2)
179  end do
180 
181 
182  ! Output is different for forward and backward simulations
183  do kz=1,numzgrid
184  do jy=0,numygrid-1
185  do ix=0,numxgrid-1
186  if (ldirect.eq.1) then
187  factor3d(ix,jy,kz)=1.e12/volume(ix,jy,kz)/outnum
188  else
189 #if defined WITH_CTBTO_PATCHES
190  factor3d(ix,jy,kz)=1.e12/volume(ix,jy,kz)/outnum
191 #else
192  factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
193 #endif
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  gridtotal=0.
205  gridsigmatotal=0.
206  gridtotalunc=0.
207  wetgridtotal=0.
208  wetgridsigmatotal=0.
209  wetgridtotalunc=0.
210  drygridtotal=0.
211  drygridsigmatotal=0.
212  drygridtotalunc=0.
213 
214  do ks=1,nspec
215 
216  write(anspec,'(i3.3)') ks
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_'//adate// &
220  atime//'_'//anspec,form='unformatted')
221  else
222  open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_'//adate// &
223  atime//'_'//anspec,form='unformatted')
224  endif
225  write(unitoutgrid) itime
226  endif
227 
228  if ((iout.eq.2).or.(iout.eq.3)) then ! mixing ratio
229  open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_'//adate// &
230  atime//'_'//anspec,form='unformatted')
231 
232  write(unitoutgridppt) itime
233  endif
234 
235  do kp=1,maxpointspec_act
236  do nage=1,nageclass
237 
238  do jy=0,numygrid-1
239  do ix=0,numxgrid-1
240 
241  ! WET DEPOSITION
242  if ((wetdep).and.(ldirect.gt.0)) then
243  do l=1,nclassunc
244  auxgrid(l)=wetgridunc(ix,jy,ks,kp,l,nage)
245  end do
246  call mean(auxgrid,wetgrid(ix,jy), &
247  wetgridsigma(ix,jy),nclassunc)
248  ! Multiply by number of classes to get total concentration
249  wetgrid(ix,jy)=wetgrid(ix,jy) &
250  *nclassunc
251  wetgridtotal=wetgridtotal+wetgrid(ix,jy)
252  ! Calculate standard deviation of the mean
253  wetgridsigma(ix,jy)= &
254  wetgridsigma(ix,jy)* &
255  sqrt(real(nclassunc))
256  wetgridsigmatotal=wetgridsigmatotal+ &
257  wetgridsigma(ix,jy)
258  endif
259 
260  ! DRY DEPOSITION
261  if ((drydep).and.(ldirect.gt.0)) then
262  do l=1,nclassunc
263  auxgrid(l)=drygridunc(ix,jy,ks,kp,l,nage)
264  end do
265  call mean(auxgrid,drygrid(ix,jy), &
266  drygridsigma(ix,jy),nclassunc)
267  ! Multiply by number of classes to get total concentration
268  drygrid(ix,jy)=drygrid(ix,jy)* &
269  nclassunc
270  drygridtotal=drygridtotal+drygrid(ix,jy)
271  ! Calculate standard deviation of the mean
272  drygridsigma(ix,jy)= &
273  drygridsigma(ix,jy)* &
274  sqrt(real(nclassunc))
275 125 drygridsigmatotal=drygridsigmatotal+ &
276  drygridsigma(ix,jy)
277  endif
278 
279  ! CONCENTRATION OR MIXING RATIO
280  do kz=1,numzgrid
281  do l=1,nclassunc
282  auxgrid(l)=gridunc(ix,jy,kz,ks,kp,l,nage)
283  end do
284  call mean(auxgrid,grid(ix,jy,kz), &
285  gridsigma(ix,jy,kz),nclassunc)
286  ! Multiply by number of classes to get total concentration
287  grid(ix,jy,kz)= &
288  grid(ix,jy,kz)*nclassunc
289  gridtotal=gridtotal+grid(ix,jy,kz)
290  ! Calculate standard deviation of the mean
291  gridsigma(ix,jy,kz)= &
292  gridsigma(ix,jy,kz)* &
293  sqrt(real(nclassunc))
294  gridsigmatotal=gridsigmatotal+ &
295  gridsigma(ix,jy,kz)
296  end do
297  end do
298  end do
299 
300 
301 
302 
303  !*******************************************************************
304  ! Generate output: may be in concentration (ng/m3) or in mixing
305  ! ratio (ppt) or both
306  ! Output the position and the values alternated multiplied by
307  ! 1 or -1, first line is number of values, number of positions
308  ! For backward simulations, the unit is seconds, stored in grid_time
309  !*******************************************************************
310 
311  ! Concentration output
312  !*********************
313  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
314 
315  ! Wet deposition
316  sp_count_i=0
317  sp_count_r=0
318  sp_fact=-1.
319  sp_zer=.true.
320  if ((ldirect.eq.1).and.(wetdep)) then
321  do jy=0,numygrid-1
322  do ix=0,numxgrid-1
323  !oncentraion greater zero
324  if (wetgrid(ix,jy).gt.smallnum) then
325  if (sp_zer.eqv..true.) then ! first non zero value
326  sp_count_i=sp_count_i+1
327  sparse_dump_i(sp_count_i)=ix+jy*numxgrid
328  sp_zer=.false.
329  sp_fact=sp_fact*(-1.)
330  endif
331  sp_count_r=sp_count_r+1
332  sparse_dump_r(sp_count_r)= &
333  sp_fact*1.e12*wetgrid(ix,jy)/area(ix,jy)
334  ! sparse_dump_u(sp_count_r)=
335  !+ 1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
336  else ! concentration is zero
337  sp_zer=.true.
338  endif
339  end do
340  end do
341  else
342  sp_count_i=0
343  sp_count_r=0
344  endif
345  write(unitoutgrid) sp_count_i
346  write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
347  write(unitoutgrid) sp_count_r
348  write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
349  ! write(unitoutgrid) sp_count_u
350  ! write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
351 
352  ! Dry deposition
353  sp_count_i=0
354  sp_count_r=0
355  sp_fact=-1.
356  sp_zer=.true.
357  if ((ldirect.eq.1).and.(drydep)) then
358  do jy=0,numygrid-1
359  do ix=0,numxgrid-1
360  if (drygrid(ix,jy).gt.smallnum) then
361  if (sp_zer.eqv..true.) then ! first non zero value
362  sp_count_i=sp_count_i+1
363  sparse_dump_i(sp_count_i)=ix+jy*numxgrid
364  sp_zer=.false.
365  sp_fact=sp_fact*(-1.)
366  endif
367  sp_count_r=sp_count_r+1
368  sparse_dump_r(sp_count_r)= &
369  sp_fact* &
370  1.e12*drygrid(ix,jy)/area(ix,jy)
371  ! sparse_dump_u(sp_count_r)=
372  !+ 1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
373  else ! concentration is zero
374  sp_zer=.true.
375  endif
376  end do
377  end do
378  else
379  sp_count_i=0
380  sp_count_r=0
381  endif
382  write(unitoutgrid) sp_count_i
383  write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
384  write(unitoutgrid) sp_count_r
385  write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
386  ! write(*,*) sp_count_u
387  ! write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
388 
389 
390 
391  ! Concentrations
392  sp_count_i=0
393  sp_count_r=0
394  sp_fact=-1.
395  sp_zer=.true.
396  do kz=1,numzgrid
397  do jy=0,numygrid-1
398  do ix=0,numxgrid-1
399  if (grid(ix,jy,kz).gt.smallnum) then
400  if (sp_zer.eqv..true.) then ! first non zero value
401  sp_count_i=sp_count_i+1
402  sparse_dump_i(sp_count_i)= &
403  ix+jy*numxgrid+kz*numxgrid*numygrid
404  sp_zer=.false.
405  sp_fact=sp_fact*(-1.)
406  endif
407  sp_count_r=sp_count_r+1
408  sparse_dump_r(sp_count_r)= &
409  sp_fact* &
410  grid(ix,jy,kz)* &
411  factor3d(ix,jy,kz)/tot_mu(ks,kp)
412  ! if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
413  ! + write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
414  ! sparse_dump_u(sp_count_r)=
415  !+ ,gridsigma(ix,jy,kz,ks,kp,nage)*
416  !+ factor(ix,jy,kz)/tot_mu(ks,kp)
417  else ! concentration is zero
418  sp_zer=.true.
419  endif
420  end do
421  end do
422  end do
423  write(unitoutgrid) sp_count_i
424  write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
425  write(unitoutgrid) sp_count_r
426  write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
427  ! write(unitoutgrid) sp_count_u
428  ! write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
429 
430 
431 
432  endif ! concentration output
433 
434  ! Mixing ratio output
435  !********************
436 
437  if ((iout.eq.2).or.(iout.eq.3)) then ! mixing ratio
438 
439  ! Wet deposition
440  sp_count_i=0
441  sp_count_r=0
442  sp_fact=-1.
443  sp_zer=.true.
444  if ((ldirect.eq.1).and.(wetdep)) then
445  do jy=0,numygrid-1
446  do ix=0,numxgrid-1
447  if (wetgrid(ix,jy).gt.smallnum) then
448  if (sp_zer.eqv..true.) then ! first non zero value
449  sp_count_i=sp_count_i+1
450  sparse_dump_i(sp_count_i)= &
451  ix+jy*numxgrid
452  sp_zer=.false.
453  sp_fact=sp_fact*(-1.)
454  endif
455  sp_count_r=sp_count_r+1
456  sparse_dump_r(sp_count_r)= &
457  sp_fact* &
458  1.e12*wetgrid(ix,jy)/area(ix,jy)
459  ! sparse_dump_u(sp_count_r)=
460  ! + ,1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
461  else ! concentration is zero
462  sp_zer=.true.
463  endif
464  end do
465  end do
466  else
467  sp_count_i=0
468  sp_count_r=0
469  endif
470  write(unitoutgridppt) sp_count_i
471  write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
472  write(unitoutgridppt) sp_count_r
473  write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
474  ! write(unitoutgridppt) sp_count_u
475  ! write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
476 
477 
478  ! Dry deposition
479  sp_count_i=0
480  sp_count_r=0
481  sp_fact=-1.
482  sp_zer=.true.
483  if ((ldirect.eq.1).and.(drydep)) then
484  do jy=0,numygrid-1
485  do ix=0,numxgrid-1
486  if (drygrid(ix,jy).gt.smallnum) then
487  if (sp_zer.eqv..true.) then ! first non zero value
488  sp_count_i=sp_count_i+1
489  sparse_dump_i(sp_count_i)= &
490  ix+jy*numxgrid
491  sp_zer=.false.
492  sp_fact=sp_fact*(-1)
493  endif
494  sp_count_r=sp_count_r+1
495  sparse_dump_r(sp_count_r)= &
496  sp_fact* &
497  1.e12*drygrid(ix,jy)/area(ix,jy)
498  ! sparse_dump_u(sp_count_r)=
499  ! + ,1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
500  else ! concentration is zero
501  sp_zer=.true.
502  endif
503  end do
504  end do
505  else
506  sp_count_i=0
507  sp_count_r=0
508  endif
509  write(unitoutgridppt) sp_count_i
510  write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
511  write(unitoutgridppt) sp_count_r
512  write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
513  ! write(unitoutgridppt) sp_count_u
514  ! write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
515 
516 
517  ! Mixing ratios
518  sp_count_i=0
519  sp_count_r=0
520  sp_fact=-1.
521  sp_zer=.true.
522  do kz=1,numzgrid
523  do jy=0,numygrid-1
524  do ix=0,numxgrid-1
525  if (grid(ix,jy,kz).gt.smallnum) then
526  if (sp_zer.eqv..true.) then ! first non zero value
527  sp_count_i=sp_count_i+1
528  sparse_dump_i(sp_count_i)= &
529  ix+jy*numxgrid+kz*numxgrid*numygrid
530  sp_zer=.false.
531  sp_fact=sp_fact*(-1.)
532  endif
533  sp_count_r=sp_count_r+1
534  sparse_dump_r(sp_count_r)= &
535  sp_fact* &
536  1.e12*grid(ix,jy,kz) &
537  /volume(ix,jy,kz)/outnum* &
538  weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
539  ! sparse_dump_u(sp_count_r)=
540  !+ ,1.e12*gridsigma(ix,jy,kz,ks,kp,nage)/volume(ix,jy,kz)/
541  !+ outnum*weightair/weightmolar(ks)/
542  !+ densityoutgrid(ix,jy,kz)
543  else ! concentration is zero
544  sp_zer=.true.
545  endif
546  end do
547  end do
548  end do
549  write(unitoutgridppt) sp_count_i
550  write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
551  write(unitoutgridppt) sp_count_r
552  write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
553  ! write(unitoutgridppt) sp_count_u
554  ! write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
555 
556  endif ! output for ppt
557 
558  end do
559  end do
560 
561  close(unitoutgridppt)
562  close(unitoutgrid)
563 
564  end do
565 
566  if (gridtotal.gt.0.) gridtotalunc=gridsigmatotal/gridtotal
567  if (wetgridtotal.gt.0.) wetgridtotalunc=wetgridsigmatotal/ &
568  wetgridtotal
569  if (drygridtotal.gt.0.) drygridtotalunc=drygridsigmatotal/ &
570  drygridtotal
571 
572  ! Dump of receptor concentrations
573 
574  if (numreceptor.gt.0 .and. (iout.eq.2 .or. iout.eq.3) ) then
575  write(unitoutreceptppt) itime
576  do ks=1,nspec
577  write(unitoutreceptppt) (1.e12*creceptor(i,ks)/outnum* &
578  weightair/weightmolar(ks)/densityoutrecept(i),i=1,numreceptor)
579  end do
580  endif
581 
582  ! Dump of receptor concentrations
583 
584  if (numreceptor.gt.0) then
585  write(unitoutrecept) itime
586  do ks=1,nspec
587  write(unitoutrecept) (1.e12*creceptor(i,ks)/outnum, &
588  i=1,numreceptor)
589  end do
590  endif
591 
592 
593 
594  ! Reinitialization of grid
595  !*************************
596 
597  do ks=1,nspec
598  do kp=1,maxpointspec_act
599  do i=1,numreceptor
600  creceptor(i,ks)=0.
601  end do
602  do jy=0,numygrid-1
603  do ix=0,numxgrid-1
604  do l=1,nclassunc
605  do nage=1,nageclass
606  do kz=1,numzgrid
607  gridunc(ix,jy,kz,ks,kp,l,nage)=0.
608  end do
609  end do
610  end do
611  end do
612  end do
613  end do
614  end do
615 
616 
617 end subroutine concoutput
subroutine caldate(juldate, yyyymmdd, hhmiss)
Definition: caldate.f90:22
subroutine concoutput(itime, outnum, gridtotalunc, wetgridtotalunc, drygridtotalunc)
Definition: concoutput.F90:22
subroutine mean(x, xm, xs, number)
Definition: mean.f90:22