source: flexpart.git/src/conccalc.f90 @ 5f2e8f6

flexpart-noresm
Last change on this file since 5f2e8f6 was 5f2e8f6, checked in by Ignacio Pisso <Ignacio.Pisso@…>, 8 years ago

new flexpart noresm code in src

  • Property mode set to 100755
File size: 14.9 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  ! For forward simulations, make a loop over the number of species;
64  ! for backward simulations, make an additional loop over the
65  ! releasepoints
66  !***************************************************************************
67
68
69  do i=1,numpart
70    if (itra1(i).ne.itime) goto 20
71
72  ! Determine age class of the particle
73    itage=abs(itra1(i)-itramem(i))
74    do nage=1,nageclass
75      if (itage.lt.lage(nage)) goto 33
76    end do
7733   continue
78
79
80  ! For special runs, interpolate the air density to the particle position
81  !************************************************************************
82  !***********************************************************************
83  !AF IND_SOURCE switches between different units for concentrations at the source
84  !Af    NOTE that in backward simulations the release of particles takes place
85  !Af    at the receptor and the sampling at the source.
86  !Af          1="mass"
87  !Af          2="mass mixing ratio"
88  !Af IND_RECEPTOR switches between different units for concentrations at the receptor
89  !Af          1="mass"
90  !Af          2="mass mixing ratio"
91
92  !Af switches for the conccalcfile:
93  !AF IND_SAMP =  0 : xmass * 1
94  !Af IND_SAMP = -1 : xmass / rho
95
96  !Af ind_samp is defined in readcommand.f
97
98    if ( ind_samp .eq. -1 ) then
99
100      ix=int(xtra1(i))
101      jy=int(ytra1(i))
102      ixp=ix+1
103      jyp=jy+1
104      ddx=xtra1(i)-real(ix)
105      ddy=ytra1(i)-real(jy)
106      rddx=1.-ddx
107      rddy=1.-ddy
108      p1=rddx*rddy
109      p2=ddx*rddy
110      p3=rddx*ddy
111      p4=ddx*ddy
112
113      do il=2,nz
114        if (height(il).gt.ztra1(i)) then
115          indz=il-1
116          indzp=il
117          goto 6
118        endif
119      end do
1206     continue
121
122      dz1=ztra1(i)-height(indz)
123      dz2=height(indzp)-ztra1(i)
124      dz=1./(dz1+dz2)
125
126  ! Take density from 2nd wind field in memory (accurate enough, no time interpolation needed)
127  !*****************************************************************************
128      do ind=indz,indzp
129        rhoprof(ind-indz+1)=p1*rho(ix ,jy ,ind,2) &
130             +p2*rho(ixp,jy ,ind,2) &
131             +p3*rho(ix ,jyp,ind,2) &
132             +p4*rho(ixp,jyp,ind,2)
133      end do
134      rhoi=(dz1*rhoprof(2)+dz2*rhoprof(1))*dz
135   elseif (ind_samp.eq.0) then
136      rhoi = 1.
137   endif
138
139
140  !****************************************************************************
141  ! 1. Evaluate grid concentrations using a uniform kernel of bandwidths dx, dy
142  !****************************************************************************
143
144
145  ! For backward simulations, look from which release point the particle comes from
146  ! For domain-filling trajectory option, npoint contains a consecutive particle
147  ! number, not the release point information. Therefore, nrelpointer is set to 1
148  ! for the domain-filling option.
149  !*****************************************************************************
150
151    if ((ioutputforeachrelease.eq.0).or.(mdomainfill.eq.1)) then
152       nrelpointer=1
153    else
154       nrelpointer=npoint(i)
155    endif
156
157    do kz=1,numzgrid                ! determine height of cell
158      if (outheight(kz).gt.ztra1(i)) goto 21
159    end do
16021   continue
161    if (kz.le.numzgrid) then           ! inside output domain
162
163
164  !********************************
165  ! Do everything for mother domain
166  !********************************
167
168      xl=(xtra1(i)*dx+xoutshift)/dxout
169      if (xl.ge.360) xl=xl-360 !added by mc to allow output grid starting at 180 with wind field origin at -177.5 ;by mc 17/6/2014
170      yl=(ytra1(i)*dy+youtshift)/dyout
171      ix=int(xl)
172      if (xl.lt.0.) ix=ix-1
173      jy=int(yl)
174      if (yl.lt.0.) jy=jy-1
175
176  ! if (i.eq.10000) write(*,*) itime,xtra1(i),ytra1(i),ztra1(i),xl,yl
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)) 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          do ks=1,nspec
191            gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
192                 gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
193                 xmass1(i,ks)/rhoi*weight
194          end do
195        endif
196
197      else                                 ! attribution via uniform kernel
198
199        ddx=xl-real(ix)                   ! distance to left cell border
200        ddy=yl-real(jy)                   ! distance to lower cell border
201        if (ddx.gt.0.5) then
202          ixp=ix+1
203          wx=1.5-ddx
204        else
205          ixp=ix-1
206          wx=0.5+ddx
207        endif
208
209        if (ddy.gt.0.5) then
210          jyp=jy+1
211          wy=1.5-ddy
212        else
213          jyp=jy-1
214          wy=0.5+ddy
215        endif
216
217
218  ! Determine mass fractions for four grid points
219  !**********************************************
220
221        if ((ix.ge.0).and.(ix.le.numxgrid-1)) then
222          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
223            w=wx*wy
224            do ks=1,nspec
225              gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
226                   gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
227                   xmass1(i,ks)/rhoi*weight*w
228            end do
229          endif
230
231          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
232            w=wx*(1.-wy)
233            do ks=1,nspec
234              gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
235                   gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
236                   xmass1(i,ks)/rhoi*weight*w
237            end do
238          endif
239        endif
240
241
242        if ((ixp.ge.0).and.(ixp.le.numxgrid-1)) then
243          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
244            w=(1.-wx)*(1.-wy)
245            do ks=1,nspec
246              gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
247                   gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
248                   xmass1(i,ks)/rhoi*weight*w
249            end do
250          endif
251
252          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
253            w=(1.-wx)*wy
254            do ks=1,nspec
255              gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
256                   gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
257                   xmass1(i,ks)/rhoi*weight*w
258            end do
259          endif
260        endif
261      endif
262
263
264
265  !************************************
266  ! Do everything for the nested domain
267  !************************************
268
269      if (nested_output.eq.1) then
270        xl=(xtra1(i)*dx+xoutshiftn)/dxoutn
271        yl=(ytra1(i)*dy+youtshiftn)/dyoutn
272        ix=int(xl)
273        if (xl.lt.0.) ix=ix-1
274        jy=int(yl)
275        if (yl.lt.0.) jy=jy-1
276
277
278  ! For particles aged less than 3 hours, attribute particle mass to grid cell
279  ! it resides in rather than use the kernel, in order to avoid its smoothing effect.
280  ! For older particles, use the uniform kernel.
281  ! If a particle is close to the domain boundary, do not use the kernel either.
282  !*****************************************************************************
283
284        if ((itage.lt.10800).or.(xl.lt.0.5).or.(yl.lt.0.5).or. &
285             (xl.gt.real(numxgridn-1)-0.5).or. &
286             (yl.gt.real(numygridn-1)-0.5)) then             ! no kernel, direct attribution to grid cell
287          if ((ix.ge.0).and.(jy.ge.0).and.(ix.le.numxgridn-1).and. &
288               (jy.le.numygridn-1)) then
289            do ks=1,nspec
290              griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
291                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
292                   xmass1(i,ks)/rhoi*weight
293            end do
294          endif
295
296        else                                 ! attribution via uniform kernel
297
298          ddx=xl-real(ix)                   ! distance to left cell border
299          ddy=yl-real(jy)                   ! distance to lower cell border
300          if (ddx.gt.0.5) then
301            ixp=ix+1
302            wx=1.5-ddx
303          else
304            ixp=ix-1
305            wx=0.5+ddx
306          endif
307
308          if (ddy.gt.0.5) then
309            jyp=jy+1
310            wy=1.5-ddy
311          else
312            jyp=jy-1
313            wy=0.5+ddy
314          endif
315
316
317  ! Determine mass fractions for four grid points
318  !**********************************************
319
320          if ((ix.ge.0).and.(ix.le.numxgridn-1)) then
321            if ((jy.ge.0).and.(jy.le.numygridn-1)) then
322              w=wx*wy
323              do ks=1,nspec
324                griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
325                     griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
326                     xmass1(i,ks)/rhoi*weight*w
327              end do
328            endif
329
330            if ((jyp.ge.0).and.(jyp.le.numygridn-1)) then
331              w=wx*(1.-wy)
332              do ks=1,nspec
333                griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
334                     griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
335                     xmass1(i,ks)/rhoi*weight*w
336              end do
337            endif
338          endif
339
340
341          if ((ixp.ge.0).and.(ixp.le.numxgridn-1)) then
342            if ((jyp.ge.0).and.(jyp.le.numygridn-1)) then
343              w=(1.-wx)*(1.-wy)
344              do ks=1,nspec
345                griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
346                     griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
347                     xmass1(i,ks)/rhoi*weight*w
348              end do
349            endif
350
351            if ((jy.ge.0).and.(jy.le.numygridn-1)) then
352              w=(1.-wx)*wy
353              do ks=1,nspec
354                griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
355                     griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
356                     xmass1(i,ks)/rhoi*weight*w
357              end do
358            endif
359          endif
360        endif
361
362      endif
363    endif
36420  continue
365  end do
366
367  !***********************************************************************
368  ! 2. Evaluate concentrations at receptor points, using the kernel method
369  !***********************************************************************
370
371  do n=1,numreceptor
372
373
374  ! Reset concentrations
375  !*********************
376
377    do ks=1,nspec
378      c(ks)=0.
379    end do
380
381
382  ! Estimate concentration at receptor
383  !***********************************
384
385    do i=1,numpart
386
387      if (itra1(i).ne.itime) goto 40
388      itage=abs(itra1(i)-itramem(i))
389
390      hz=min(50.+0.3*sqrt(real(itage)),hzmax)
391      zd=ztra1(i)/hz
392      if (zd.gt.1.) goto 40          ! save computing time, leave loop
393
394      hx=min((0.29+2.222e-3*sqrt(real(itage)))*dx+ &
395           real(itage)*1.2e-5,hxmax)                     ! 80 km/day
396      xd=(xtra1(i)-xreceptor(n))/hx
397      if (xd*xd.gt.1.) goto 40       ! save computing time, leave loop
398
399      hy=min((0.18+1.389e-3*sqrt(real(itage)))*dy+ &
400           real(itage)*7.5e-6,hymax)                     ! 80 km/day
401      yd=(ytra1(i)-yreceptor(n))/hy
402      if (yd*yd.gt.1.) goto 40       ! save computing time, leave loop
403      h=hx*hy*hz
404
405      r2=xd*xd+yd*yd+zd*zd
406      if (r2.lt.1.) then
407        xkern=factor*(1.-r2)
408        do ks=1,nspec
409          c(ks)=c(ks)+xmass1(i,ks)*xkern/h
410        end do
411      endif
41240    continue
413    end do
414
415    do ks=1,nspec
416      creceptor(n,ks)=creceptor(n,ks)+2.*weight*c(ks)/receptorarea(n)
417    end do
418  end do
419
420end subroutine conccalc
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG