source: flexpart.git/src/concoutput_nest.f90 @ 4c64400

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

Completed handling of nested wind fields with cloud water (for wet deposition).

  • Property mode set to 100644
File size: 20.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  !*****************************************************************************
49  !                                                                            *
50  ! Variables:                                                                 *
51  ! outnum          number of samples                                          *
52  ! ncells          number of cells with non-zero concentrations               *
53  ! sparse          .true. if in sparse matrix format, else .false.            *
54  ! tot_mu          1 for forward, initial mass mixing ration for backw. runs  *
55  !                                                                            *
56  !*****************************************************************************
57
58  use unc_mod
59  use point_mod
60  use outg_mod
61  use par_mod
62  use com_mod
63  use mean_mod
64
65  implicit none
66
67  real(kind=dp) :: jul
68  integer :: itime,i,ix,jy,kz,ks,kp,l,iix,jjy,kzz,nage,jjjjmmdd,ihmmss
69  integer :: sp_count_i,sp_count_r
70  real :: sp_fact
71  real :: outnum,densityoutrecept(maxreceptor),xl,yl
72
73  !real densityoutgrid(0:numxgrid-1,0:numygrid-1,numzgrid),
74  !    +grid(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,maxpointspec_act,
75  !    +    maxageclass)
76  !real wetgrid(0:numxgrid-1,0:numygrid-1,maxspec,maxpointspec_act,
77  !    +       maxageclass)
78  !real drygrid(0:numxgrid-1,0:numygrid-1,maxspec,
79  !    +       maxpointspec_act,maxageclass)
80  !real gridsigma(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,
81  !    +       maxpointspec_act,maxageclass),
82  !    +     drygridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
83  !    +     maxpointspec_act,maxageclass),
84  !    +     wetgridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
85  !    +     maxpointspec_act,maxageclass)
86  !real factor(0:numxgrid-1,0:numygrid-1,numzgrid)
87  !real sparse_dump_r(numxgrid*numygrid*numzgrid)
88  !integer sparse_dump_i(numxgrid*numygrid*numzgrid)
89
90  !real sparse_dump_u(numxgrid*numygrid*numzgrid)
91  real(dep_prec) :: auxgrid(nclassunc)
92  real :: halfheight,dz,dz1,dz2,tot_mu(maxspec,maxpointspec_act)
93  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
94  real,parameter :: weightair=28.97
95  logical :: sp_zer
96  character :: adate*8,atime*6
97  character(len=3) :: anspec
98  integer :: mind
99! mind        eso:added to ensure identical results between 2&3-fields versions
100
101
102  ! Determine current calendar date, needed for the file name
103  !**********************************************************
104
105  jul=bdate+real(itime,kind=dp)/86400._dp
106  call caldate(jul,jjjjmmdd,ihmmss)
107  write(adate,'(i8.8)') jjjjmmdd
108  write(atime,'(i6.6)') ihmmss
109
110
111  ! For forward simulations, output fields have dimension MAXSPEC,
112  ! for backward simulations, output fields have dimension MAXPOINT.
113  ! Thus, make loops either about nspec, or about numpoint
114  !*****************************************************************
115
116
117    if (ldirect.eq.1) then
118       do ks=1,nspec
119         do kp=1,maxpointspec_act
120           tot_mu(ks,kp)=1
121         end do
122       end do
123   else
124      do ks=1,nspec
125             do kp=1,maxpointspec_act
126               tot_mu(ks,kp)=xmass(kp,ks)
127             end do
128      end do
129    endif
130
131
132  !*******************************************************************
133  ! Compute air density: sufficiently accurate to take it
134  ! from coarse grid at some time
135  ! Determine center altitude of output layer, and interpolate density
136  ! data to that altitude
137  !*******************************************************************
138
139  mind=memind(2)
140  do kz=1,numzgrid
141    if (kz.eq.1) then
142      halfheight=outheight(1)/2.
143    else
144      halfheight=(outheight(kz)+outheight(kz-1))/2.
145    endif
146    do kzz=2,nz
147      if ((height(kzz-1).lt.halfheight).and. &
148           (height(kzz).gt.halfheight)) goto 46
149    end do
15046   kzz=max(min(kzz,nz),2)
151    dz1=halfheight-height(kzz-1)
152    dz2=height(kzz)-halfheight
153    dz=dz1+dz2
154    do jy=0,numygridn-1
155      do ix=0,numxgridn-1
156        xl=outlon0n+real(ix)*dxoutn
157        yl=outlat0n+real(jy)*dyoutn
158        xl=(xl-xlon0)/dx
159        yl=(yl-ylat0)/dy
160        iix=max(min(nint(xl),nxmin1),0)
161        jjy=max(min(nint(yl),nymin1),0)
162        ! densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
163        !      rho(iix,jjy,kzz-1,2)*dz2)/dz
164        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,mind)*dz1+ &
165             rho(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  end do
178
179
180  ! Output is different for forward and backward simulations
181    do kz=1,numzgrid
182      do jy=0,numygridn-1
183        do ix=0,numxgridn-1
184          if (ldirect.eq.1) then
185            factor3d(ix,jy,kz)=1.e12/volumen(ix,jy,kz)/outnum
186          else
187            factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
188          endif
189        end do
190      end do
191    end do
192
193  !*********************************************************************
194  ! Determine the standard deviation of the mean concentration or mixing
195  ! ratio (uncertainty of the output) and the dry and wet deposition
196  !*********************************************************************
197
198  do ks=1,nspec
199
200  write(anspec,'(i3.3)') ks
201  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
202    if (ldirect.eq.1) then
203      open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_nest_' &
204           //adate// &
205           atime//'_'//anspec,form='unformatted')
206    else
207      open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_nest_' &
208           //adate// &
209           atime//'_'//anspec,form='unformatted')
210    endif
211     write(unitoutgrid) itime
212   endif
213
214  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
215   open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_nest_' &
216        //adate// &
217        atime//'_'//anspec,form='unformatted')
218
219    write(unitoutgridppt) itime
220  endif
221
222  do kp=1,maxpointspec_act
223  do nage=1,nageclass
224
225    do jy=0,numygridn-1
226      do ix=0,numxgridn-1
227
228  ! WET DEPOSITION
229        if ((WETDEP).and.(ldirect.gt.0)) then
230            do l=1,nclassunc
231              auxgrid(l)=wetgriduncn(ix,jy,ks,kp,l,nage)
232            end do
233            call mean(auxgrid,wetgrid(ix,jy), &
234                 wetgridsigma(ix,jy),nclassunc)
235  ! Multiply by number of classes to get total concentration
236            wetgrid(ix,jy)=wetgrid(ix,jy) &
237                 *nclassunc
238  ! Calculate standard deviation of the mean
239            wetgridsigma(ix,jy)= &
240                 wetgridsigma(ix,jy)* &
241                 sqrt(real(nclassunc))
242        endif
243
244  ! DRY DEPOSITION
245        if ((DRYDEP).and.(ldirect.gt.0)) then
246            do l=1,nclassunc
247              auxgrid(l)=drygriduncn(ix,jy,ks,kp,l,nage)
248            end do
249            call mean(auxgrid,drygrid(ix,jy), &
250                 drygridsigma(ix,jy),nclassunc)
251  ! Multiply by number of classes to get total concentration
252            drygrid(ix,jy)=drygrid(ix,jy)* &
253                 nclassunc
254  ! Calculate standard deviation of the mean
255            drygridsigma(ix,jy)= &
256                 drygridsigma(ix,jy)* &
257                 sqrt(real(nclassunc))
258        endif
259
260  ! CONCENTRATION OR MIXING RATIO
261        do kz=1,numzgrid
262            do l=1,nclassunc
263              auxgrid(l)=griduncn(ix,jy,kz,ks,kp,l,nage)
264            end do
265            call mean(auxgrid,grid(ix,jy,kz), &
266                 gridsigma(ix,jy,kz),nclassunc)
267  ! Multiply by number of classes to get total concentration
268            grid(ix,jy,kz)= &
269                 grid(ix,jy,kz)*nclassunc
270  ! Calculate standard deviation of the mean
271            gridsigma(ix,jy,kz)= &
272                 gridsigma(ix,jy,kz)* &
273                 sqrt(real(nclassunc))
274        end do
275      end do
276    end do
277
278
279  !*******************************************************************
280  ! Generate output: may be in concentration (ng/m3) or in mixing
281  ! ratio (ppt) or both
282  ! Output the position and the values alternated multiplied by
283  ! 1 or -1, first line is number of values, number of positions
284  ! For backward simulations, the unit is seconds, stored in grid_time
285  !*******************************************************************
286
287  ! Concentration output
288  !*********************
289  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
290
291  ! Wet deposition
292         sp_count_i=0
293         sp_count_r=0
294         sp_fact=-1.
295         sp_zer=.true.
296         if ((ldirect.eq.1).and.(WETDEP)) then
297         do jy=0,numygridn-1
298            do ix=0,numxgridn-1
299  !oncentraion greater zero
300              if (wetgrid(ix,jy).gt.smallnum) then
301                 if (sp_zer.eqv..true.) then ! first non zero value
302                    sp_count_i=sp_count_i+1
303                    sparse_dump_i(sp_count_i)=ix+jy*numxgridn
304                    sp_zer=.false.
305                    sp_fact=sp_fact*(-1.)
306                 endif
307                 sp_count_r=sp_count_r+1
308                 sparse_dump_r(sp_count_r)= &
309                      sp_fact*1.e12*wetgrid(ix,jy)/arean(ix,jy)
310  !                sparse_dump_u(sp_count_r)=
311  !+                1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
312              else ! concentration is zero
313                  sp_zer=.true.
314              endif
315            end do
316         end do
317         else
318            sp_count_i=0
319            sp_count_r=0
320         endif
321         write(unitoutgrid) sp_count_i
322         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
323         write(unitoutgrid) sp_count_r
324         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
325  !       write(unitoutgrid) sp_count_u
326  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
327
328  ! Dry deposition
329         sp_count_i=0
330         sp_count_r=0
331         sp_fact=-1.
332         sp_zer=.true.
333         if ((ldirect.eq.1).and.(DRYDEP)) then
334          do jy=0,numygridn-1
335            do ix=0,numxgridn-1
336              if (drygrid(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* &
346                      1.e12*drygrid(ix,jy)/arean(ix,jy)
347  !                sparse_dump_u(sp_count_r)=
348  !+                1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
349              else ! concentration is zero
350                  sp_zer=.true.
351              endif
352            end do
353          end do
354         else
355            sp_count_i=0
356            sp_count_r=0
357         endif
358         write(unitoutgrid) sp_count_i
359         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
360         write(unitoutgrid) sp_count_r
361         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
362  !       write(*,*) sp_count_u
363  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
364
365
366
367  ! Concentrations
368         sp_count_i=0
369         sp_count_r=0
370         sp_fact=-1.
371         sp_zer=.true.
372          do kz=1,numzgrid
373            do jy=0,numygridn-1
374              do ix=0,numxgridn-1
375                if (grid(ix,jy,kz).gt.smallnum) then
376                  if (sp_zer.eqv..true.) then ! first non zero value
377                    sp_count_i=sp_count_i+1
378                    sparse_dump_i(sp_count_i)= &
379                         ix+jy*numxgridn+kz*numxgridn*numygridn
380                    sp_zer=.false.
381                    sp_fact=sp_fact*(-1.)
382                   endif
383                   sp_count_r=sp_count_r+1
384                   sparse_dump_r(sp_count_r)= &
385                        sp_fact* &
386                        grid(ix,jy,kz)* &
387                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
388  !                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
389  !    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
390  !                sparse_dump_u(sp_count_r)=
391  !+               ,gridsigma(ix,jy,kz,ks,kp,nage)*
392  !+               factor(ix,jy,kz)/tot_mu(ks,kp)
393              else ! concentration is zero
394                  sp_zer=.true.
395              endif
396              end do
397            end do
398          end do
399         write(unitoutgrid) sp_count_i
400         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
401         write(unitoutgrid) sp_count_r
402         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
403  !       write(unitoutgrid) sp_count_u
404  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
405
406
407
408    endif !  concentration output
409
410  ! Mixing ratio output
411  !********************
412
413  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
414
415  ! Wet deposition
416         sp_count_i=0
417         sp_count_r=0
418         sp_fact=-1.
419         sp_zer=.true.
420         if ((ldirect.eq.1).and.(WETDEP)) then
421          do jy=0,numygridn-1
422            do ix=0,numxgridn-1
423                if (wetgrid(ix,jy).gt.smallnum) then
424                  if (sp_zer.eqv..true.) then ! first non zero value
425                    sp_count_i=sp_count_i+1
426                    sparse_dump_i(sp_count_i)= &
427                         ix+jy*numxgridn
428                    sp_zer=.false.
429                    sp_fact=sp_fact*(-1.)
430                 endif
431                 sp_count_r=sp_count_r+1
432                 sparse_dump_r(sp_count_r)= &
433                      sp_fact* &
434                      1.e12*wetgrid(ix,jy)/arean(ix,jy)
435  !                sparse_dump_u(sp_count_r)=
436  !    +            ,1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
437              else ! concentration is zero
438                  sp_zer=.true.
439              endif
440            end do
441          end do
442         else
443           sp_count_i=0
444           sp_count_r=0
445         endif
446         write(unitoutgridppt) sp_count_i
447         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
448         write(unitoutgridppt) sp_count_r
449         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
450  !       write(unitoutgridppt) sp_count_u
451  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
452
453
454  ! Dry deposition
455         sp_count_i=0
456         sp_count_r=0
457         sp_fact=-1.
458         sp_zer=.true.
459         if ((ldirect.eq.1).and.(DRYDEP)) then
460          do jy=0,numygridn-1
461            do ix=0,numxgridn-1
462                if (drygrid(ix,jy).gt.smallnum) then
463                  if (sp_zer.eqv..true.) then ! first non zero value
464                    sp_count_i=sp_count_i+1
465                    sparse_dump_i(sp_count_i)= &
466                         ix+jy*numxgridn
467                    sp_zer=.false.
468                    sp_fact=sp_fact*(-1)
469                 endif
470                 sp_count_r=sp_count_r+1
471                 sparse_dump_r(sp_count_r)= &
472                      sp_fact* &
473                      1.e12*drygrid(ix,jy)/arean(ix,jy)
474  !                sparse_dump_u(sp_count_r)=
475  !    +            ,1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
476              else ! concentration is zero
477                  sp_zer=.true.
478              endif
479            end do
480          end do
481         else
482           sp_count_i=0
483           sp_count_r=0
484         endif
485         write(unitoutgridppt) sp_count_i
486         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
487         write(unitoutgridppt) sp_count_r
488         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
489  !       write(unitoutgridppt) sp_count_u
490  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
491
492
493  ! Mixing ratios
494         sp_count_i=0
495         sp_count_r=0
496         sp_fact=-1.
497         sp_zer=.true.
498          do kz=1,numzgrid
499            do jy=0,numygridn-1
500              do ix=0,numxgridn-1
501                if (grid(ix,jy,kz).gt.smallnum) then
502                  if (sp_zer.eqv..true.) then ! first non zero value
503                    sp_count_i=sp_count_i+1
504                    sparse_dump_i(sp_count_i)= &
505                         ix+jy*numxgridn+kz*numxgridn*numygridn
506                    sp_zer=.false.
507                    sp_fact=sp_fact*(-1.)
508                 endif
509                 sp_count_r=sp_count_r+1
510                 sparse_dump_r(sp_count_r)= &
511                      sp_fact* &
512                      1.e12*grid(ix,jy,kz) &
513                      /volumen(ix,jy,kz)/outnum* &
514                      weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
515  !                sparse_dump_u(sp_count_r)=
516  !+              ,1.e12*gridsigma(ix,jy,kz,ks,kp,nage)/volume(ix,jy,kz)/
517  !+              outnum*weightair/weightmolar(ks)/
518  !+              densityoutgrid(ix,jy,kz)
519              else ! concentration is zero
520                  sp_zer=.true.
521              endif
522              end do
523            end do
524          end do
525         write(unitoutgridppt) sp_count_i
526         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
527         write(unitoutgridppt) sp_count_r
528         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
529  !       write(unitoutgridppt) sp_count_u
530  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
531
532      endif ! output for ppt
533
534  end do
535  end do
536
537    close(unitoutgridppt)
538    close(unitoutgrid)
539
540  end do
541
542
543
544  ! Reinitialization of grid
545  !*************************
546
547  ! do ks=1,nspec
548  ! do kp=1,maxpointspec_act
549  !   do i=1,numreceptor
550  !     creceptor(i,ks)=0.
551  !   end do
552  !   do jy=0,numygridn-1
553  !     do ix=0,numxgridn-1
554  !       do l=1,nclassunc
555  !         do nage=1,nageclass
556  !           do kz=1,numzgrid
557  !             griduncn(ix,jy,kz,ks,kp,l,nage)=0.
558  !           end do
559  !         end do
560  !       end do
561  !     end do
562  !   end do
563  ! end do
564  ! end do
565  creceptor(:,:)=0.
566  griduncn(:,:,:,:,:,:,:)=0.
567
568
569end subroutine concoutput_nest
570
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG