source: flexpart.git/src/conccalc.f90 @ 88929bf

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 88929bf was 88929bf, checked in by Sabine <sabine.eckhardt@…>, 7 years ago

conccalc - only write pos values of xscav_frac1

  • Property mode set to 100644
File size: 18.3 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 conccalc(itime,weight)
23  !                      i     i
24  !*****************************************************************************
25  !                                                                            *
26  !     Calculation of the concentrations on a regular grid using volume       *
27  !     sampling                                                               *
28  !                                                                            *
29  !     Author: A. Stohl                                                       *
30  !                                                                            *
31  !     24 May 1996                                                            *
32  !                                                                            *
33  !     April 2000: Update to calculate age spectra                            *
34  !                 Bug fix to avoid negative conc. at the domain boundaries,  *
35  !                 as suggested by Petra Seibert                              *
36  !                                                                            *
37  !     2 July 2002: re-order if-statements in order to optimize CPU time      *
38  !                                                                            *
39  !                                                                            *
40  !*****************************************************************************
41  !                                                                            *
42  ! Variables:                                                                 *
43  ! nspeciesdim     = nspec for forward runs, 1 for backward runs              *
44  !                                                                            *
45  !*****************************************************************************
46
47  use unc_mod
48  use outg_mod
49  use par_mod
50  use com_mod
51
52  implicit none
53
54  integer :: itime,itage,i,ix,jy,ixp,jyp,kz,ks,n,nage
55  integer :: il,ind,indz,indzp,nrelpointer
56  real :: rddx,rddy,p1,p2,p3,p4,dz1,dz2,dz
57  real :: weight,hx,hy,hz,h,xd,yd,zd,xkern,r2,c(maxspec),ddx,ddy
58  real :: rhoprof(2),rhoi
59  real :: xl,yl,wx,wy,w
60  real,parameter :: factor=.596831, hxmax=6.0, hymax=4.0, hzmax=150.
61
62
63  if (.not.usekernel) then
64     write (*,*) 'NOT USING THE KERNEL!'
65  endif
66  ! For forward simulations, make a loop over the number of species;
67  ! for backward simulations, make an additional loop over the
68  ! releasepoints
69  !***************************************************************************
70
71  do i=1,numpart
72    if (itra1(i).ne.itime) goto 20
73
74  ! Determine age class of the particle
75    itage=abs(itra1(i)-itramem(i))
76    do nage=1,nageclass
77      if (itage.lt.lage(nage)) goto 33
78    end do
7933   continue
80
81
82  ! For special runs, interpolate the air density to the particle position
83  !************************************************************************
84  !***********************************************************************
85  !AF IND_SOURCE switches between different units for concentrations at the source
86  !Af    NOTE that in backward simulations the release of particles takes place
87  !Af    at the receptor and the sampling at the source.
88  !Af          1="mass"
89  !Af          2="mass mixing ratio"
90  !Af IND_RECEPTOR switches between different units for concentrations at the receptor
91  !Af          1="mass"
92  !Af          2="mass mixing ratio"
93
94  !Af switches for the conccalcfile:
95  !AF IND_SAMP =  0 : xmass * 1
96  !Af IND_SAMP = -1 : xmass / rho
97
98  !Af ind_samp is defined in readcommand.f
99
100    if ( ind_samp .eq. -1 ) then
101
102      ix=int(xtra1(i))
103      jy=int(ytra1(i))
104      ixp=ix+1
105      jyp=jy+1
106      ddx=xtra1(i)-real(ix)
107      ddy=ytra1(i)-real(jy)
108      rddx=1.-ddx
109      rddy=1.-ddy
110      p1=rddx*rddy
111      p2=ddx*rddy
112      p3=rddx*ddy
113      p4=ddx*ddy
114
115      do il=2,nz
116        if (height(il).gt.ztra1(i)) then
117          indz=il-1
118          indzp=il
119          goto 6
120        endif
121      end do
1226     continue
123
124      dz1=ztra1(i)-height(indz)
125      dz2=height(indzp)-ztra1(i)
126      dz=1./(dz1+dz2)
127
128  ! Take density from 2nd wind field in memory (accurate enough, no time interpolation needed)
129  !*****************************************************************************
130      do ind=indz,indzp
131        rhoprof(ind-indz+1)=p1*rho(ix ,jy ,ind,2) &
132             +p2*rho(ixp,jy ,ind,2) &
133             +p3*rho(ix ,jyp,ind,2) &
134             +p4*rho(ixp,jyp,ind,2)
135      end do
136      rhoi=(dz1*rhoprof(2)+dz2*rhoprof(1))*dz
137   elseif (ind_samp.eq.0) then
138      rhoi = 1.
139   endif
140
141
142  !****************************************************************************
143  ! 1. Evaluate grid concentrations using a uniform kernel of bandwidths dx, dy
144  !****************************************************************************
145
146
147  ! For backward simulations, look from which release point the particle comes from
148  ! For domain-filling trajectory option, npoint contains a consecutive particle
149  ! number, not the release point information. Therefore, nrelpointer is set to 1
150  ! for the domain-filling option.
151  !*****************************************************************************
152
153    if ((ioutputforeachrelease.eq.0).or.(mdomainfill.eq.1)) then
154       nrelpointer=1
155    else
156       nrelpointer=npoint(i)
157    endif
158
159    do kz=1,numzgrid                ! determine height of cell
160      if (outheight(kz).gt.ztra1(i)) goto 21
161    end do
16221   continue
163    if (kz.le.numzgrid) then           ! inside output domain
164
165
166  !********************************
167  ! Do everything for mother domain
168  !********************************
169
170      xl=(xtra1(i)*dx+xoutshift)/dxout
171      yl=(ytra1(i)*dy+youtshift)/dyout
172      ix=int(xl)
173      if (xl.lt.0.) ix=ix-1
174      jy=int(yl)
175      if (yl.lt.0.) jy=jy-1
176
177
178
179  ! For particles aged less than 3 hours, attribute particle mass to grid cell
180  ! it resides in rather than use the kernel, in order to avoid its smoothing effect.
181  ! For older particles, use the uniform kernel.
182  ! If a particle is close to the domain boundary, do not use the kernel either.
183  !*****************************************************************************
184
185      if (((itage.lt.10800).or.(xl.lt.0.5).or.(yl.lt.0.5).or. &
186           (xl.gt.real(numxgrid-1)-0.5).or. &
187           (yl.gt.real(numygrid-1)-0.5)).or.(.not.usekernel)) then             ! no kernel, direct attribution to grid cell
188        if ((ix.ge.0).and.(jy.ge.0).and.(ix.le.numxgrid-1).and. &
189             (jy.le.numygrid-1)) then
190          if (DRYBKDEP.or.WETBKDEP) then
191             do ks=1,nspec
192               gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
193                 gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
194                 xmass1(i,ks)/rhoi*weight*max(xscav_frac1(i,ks),0.0)
195             end do
196          else
197             do ks=1,nspec
198               gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
199                 gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
200                 xmass1(i,ks)/rhoi*weight
201             end do
202          endif
203        endif
204
205      else                                 ! attribution via uniform kernel
206
207        ddx=xl-real(ix)                   ! distance to left cell border
208        ddy=yl-real(jy)                   ! distance to lower cell border
209        if (ddx.gt.0.5) then
210          ixp=ix+1
211          wx=1.5-ddx
212        else
213          ixp=ix-1
214          wx=0.5+ddx
215        endif
216
217        if (ddy.gt.0.5) then
218          jyp=jy+1
219          wy=1.5-ddy
220        else
221          jyp=jy-1
222          wy=0.5+ddy
223        endif
224
225
226  ! Determine mass fractions for four grid points
227  !**********************************************
228
229        if ((ix.ge.0).and.(ix.le.numxgrid-1)) then
230          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
231            w=wx*wy
232            if (DRYBKDEP.or.WETBKDEP) then
233               do ks=1,nspec
234                 gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
235                   gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
236                   xmass1(i,ks)/rhoi*w*weight*max(xscav_frac1(i,ks),0.0)
237               end do
238            else
239               do ks=1,nspec
240                 gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
241                   gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
242                   xmass1(i,ks)/rhoi*weight*w
243               end do
244            endif
245          endif
246
247          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
248            w=wx*(1.-wy)
249            if (DRYBKDEP.or.WETBKDEP) then
250              do ks=1,nspec
251                 gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
252                   gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
253                   xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
254               end do
255             else
256              do ks=1,nspec
257                 gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
258                   gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
259                   xmass1(i,ks)/rhoi*weight*w
260               end do
261             endif
262          endif
263        endif !ix ge 0
264
265
266        if ((ixp.ge.0).and.(ixp.le.numxgrid-1)) then
267          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
268            w=(1.-wx)*(1.-wy)
269            if (DRYBKDEP.or.WETBKDEP) then
270               do ks=1,nspec
271                 gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
272                   gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
273                   xmass1(i,ks)/rhoi*w*weight*max(xscav_frac1(i,ks),0.0)
274               end do
275            else
276               do ks=1,nspec
277                 gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
278                   gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
279                   xmass1(i,ks)/rhoi*weight*w
280               end do
281            endif
282          endif
283
284          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
285            w=(1.-wx)*wy
286            if (DRYBKDEP.or.WETBKDEP) then
287               do ks=1,nspec
288                 gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
289                   gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
290                   xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
291               end do
292            else
293               do ks=1,nspec
294                 gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
295                   gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
296                   xmass1(i,ks)/rhoi*weight*w
297               end do
298            endif
299          endif
300        endif !ixp ge 0
301     endif
302
303  !************************************
304  ! Do everything for the nested domain
305  !************************************
306
307      if (nested_output.eq.1) then
308        xl=(xtra1(i)*dx+xoutshiftn)/dxoutn
309        yl=(ytra1(i)*dy+youtshiftn)/dyoutn
310        ix=int(xl)
311        if (xl.lt.0.) ix=ix-1
312        jy=int(yl)
313        if (yl.lt.0.) jy=jy-1
314
315
316  ! For particles aged less than 3 hours, attribute particle mass to grid cell
317  ! it resides in rather than use the kernel, in order to avoid its smoothing effect.
318  ! For older particles, use the uniform kernel.
319  ! If a particle is close to the domain boundary, do not use the kernel either.
320  !*****************************************************************************
321
322        if ((itage.lt.10800).or.(xl.lt.0.5).or.(yl.lt.0.5).or. &
323             (xl.gt.real(numxgridn-1)-0.5).or. &
324             (yl.gt.real(numygridn-1)-0.5).or.(.not.usekernel)) then             ! no kernel, direct attribution to grid cell
325          if ((ix.ge.0).and.(jy.ge.0).and.(ix.le.numxgridn-1).and. &
326               (jy.le.numygridn-1)) then
327            if (DRYBKDEP.or.WETBKDEP) then
328               do ks=1,nspec
329                 griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
330                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
331                   xmass1(i,ks)/rhoi*weight*max(xscav_frac1(i,ks),0.0)
332               end do
333            else
334               do ks=1,nspec
335                 griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
336                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
337                   xmass1(i,ks)/rhoi*weight
338               end do
339            endif
340          endif
341
342        else                                 ! attribution via uniform kernel
343
344          ddx=xl-real(ix)                   ! distance to left cell border
345          ddy=yl-real(jy)                   ! distance to lower cell border
346          if (ddx.gt.0.5) then
347            ixp=ix+1
348            wx=1.5-ddx
349          else
350            ixp=ix-1
351            wx=0.5+ddx
352          endif
353
354          if (ddy.gt.0.5) then
355            jyp=jy+1
356            wy=1.5-ddy
357          else
358            jyp=jy-1
359            wy=0.5+ddy
360          endif
361
362
363  ! Determine mass fractions for four grid points
364  !**********************************************
365
366          if ((ix.ge.0).and.(ix.le.numxgridn-1)) then
367            if ((jy.ge.0).and.(jy.le.numygridn-1)) then
368              w=wx*wy
369              if (DRYBKDEP.or.WETBKDEP) then
370                 do ks=1,nspec
371                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
372                     griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
373                     xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
374                 end do
375              else
376                do ks=1,nspec
377                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
378                     griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
379                     xmass1(i,ks)/rhoi*weight*w
380                 end do
381              endif
382            endif
383
384            if ((jyp.ge.0).and.(jyp.le.numygridn-1)) then
385              w=wx*(1.-wy)
386              if (DRYBKDEP.or.WETBKDEP) then
387                 do ks=1,nspec
388                   griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
389                     griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
390                     xmass1(i,ks)/rhoi*weight*w*xscav_frac1(i,ks)
391                 end do
392              else
393                 do ks=1,nspec
394                   griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
395                     griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
396                     xmass1(i,ks)/rhoi*weight*w
397                 end do
398              endif
399            endif
400          endif
401
402
403          if ((ixp.ge.0).and.(ixp.le.numxgridn-1)) then
404            if ((jyp.ge.0).and.(jyp.le.numygridn-1)) then
405              w=(1.-wx)*(1.-wy)
406              if (DRYBKDEP.or.WETBKDEP) then
407                 do ks=1,nspec
408                   griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
409                     griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
410                     xmass1(i,ks)/rhoi*weight*w*xscav_frac1(i,ks)
411                 end do
412              else
413                 do ks=1,nspec
414                   griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
415                     griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
416                     xmass1(i,ks)/rhoi*weight*w
417                 end do
418              endif
419            endif
420
421            if ((jy.ge.0).and.(jy.le.numygridn-1)) then
422              w=(1.-wx)*wy
423              if (DRYBKDEP.or.WETBKDEP) then
424                 do ks=1,nspec
425                   griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
426                     griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
427                     xmass1(i,ks)/rhoi*weight*w*xscav_frac1(i,ks)
428                 end do
429              else
430                 do ks=1,nspec
431                    griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
432                     griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
433                     xmass1(i,ks)/rhoi*weight*w
434                 end do
435              endif
436            endif
437          endif
438        endif
439      endif
440    endif
44120  continue
442  end do
443
444  !***********************************************************************
445  ! 2. Evaluate concentrations at receptor points, using the kernel method
446  !***********************************************************************
447
448  do n=1,numreceptor
449
450
451  ! Reset concentrations
452  !*********************
453
454    do ks=1,nspec
455      c(ks)=0.
456    end do
457
458
459  ! Estimate concentration at receptor
460  !***********************************
461
462    do i=1,numpart
463
464      if (itra1(i).ne.itime) goto 40
465      itage=abs(itra1(i)-itramem(i))
466
467      hz=min(50.+0.3*sqrt(real(itage)),hzmax)
468      zd=ztra1(i)/hz
469      if (zd.gt.1.) goto 40          ! save computing time, leave loop
470
471      hx=min((0.29+2.222e-3*sqrt(real(itage)))*dx+ &
472           real(itage)*1.2e-5,hxmax)                     ! 80 km/day
473      xd=(xtra1(i)-xreceptor(n))/hx
474      if (xd*xd.gt.1.) goto 40       ! save computing time, leave loop
475
476      hy=min((0.18+1.389e-3*sqrt(real(itage)))*dy+ &
477           real(itage)*7.5e-6,hymax)                     ! 80 km/day
478      yd=(ytra1(i)-yreceptor(n))/hy
479      if (yd*yd.gt.1.) goto 40       ! save computing time, leave loop
480      h=hx*hy*hz
481
482      r2=xd*xd+yd*yd+zd*zd
483      if (r2.lt.1.) then
484        xkern=factor*(1.-r2)
485        do ks=1,nspec
486          c(ks)=c(ks)+xmass1(i,ks)*xkern/h
487        end do
488      endif
48940    continue
490    end do
491
492    do ks=1,nspec
493      creceptor(n,ks)=creceptor(n,ks)+2.*weight*c(ks)/receptorarea(n)
494    end do
495  end do
496
497end subroutine conccalc
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG