Ticket #93: conccalc_irreg.f90

File conccalc_irreg.f90, 15.2 KB (added by jbrioude, 10 years ago)
Line 
1!***********************************************************************
2!* Copyright 2012,2013                                                *
3!* Jerome Brioude, Delia Arnold, Andreas Stohl, Wayne Angevine,       *
4!* John Burkhart, Massimo Cassiani, Adam Dingwell, Richard C Easter, Sabine Eckhardt,*
5!* Stephanie Evan, Jerome D Fast, Don Morton, Ignacio Pisso,          *
6!* Petra Seibert, Gerard Wotawa, Caroline Forster, Harald Sodemann,   *
7!*                                                                     *
8!* This file is part of FLEXPART WRF                                   *
9!*                                                                     *
10!* FLEXPART is free software: you can redistribute it and/or modify    *
11!* it under the terms of the GNU General Public License as published by*
12!* the Free Software Foundation, either version 3 of the License, or   *
13!* (at your option) any later version.                                 *
14!*                                                                     *
15!* FLEXPART is distributed in the hope that it will be useful,         *
16!* but WITHOUT ANY WARRANTY; without even the implied warranty of      *
17!* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
18!* GNU General Public License for more details.                        *
19!*                                                                     *
20!* You should have received a copy of the GNU General Public License   *
21!* along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
22!***********************************************************************
23
24subroutine conccalc_irreg(itime,weight)
25  !                      i     i
26  !*****************************************************************************
27  !                                                                            *
28  !     Calculation of the concentrations on a regular grid using volume       *
29  !     sampling                                                               *
30  !                                                                            *
31  !     Author: A. Stohl                                                       *
32  !                                                                            *
33  !     24 May 1996                                                            *
34  !                                                                            *
35  !     April 2000: Update to calculate age spectra                            *
36  !                 Bug fix to avoid negative conc. at the domain boundaries,  *
37  !                 as suggested by Petra Seibert                              *
38  !                                                                            *
39  !     2 July 2002: re-order if-statements in order to optimize CPU time      *
40  !                                                                            *
41  !                                                                            *
42  !*****************************************************************************
43  !                                                                            *
44  ! Variables:                                                                 *
45  ! nspeciesdim     = nspec for forward runs, 1 for backward runs              *
46  !                                                                            *
47  !*****************************************************************************
48
49  use unc_mod
50  use outg_mod
51  use par_mod
52  use com_mod
53
54  implicit none
55
56  integer :: itime,itage,i,ix,jy,ixp,jyp,kz,ks,n,nage
57  integer :: il,ind,indz,indzp,nrelpointer,ii,temp,incl
58  real :: rddx,rddy,p1,p2,p3,p4,dz1,dz2,dz
59  real :: weight,hx,hy,hz,h,xd,yd,zd,xkern,r2,c(maxspec),ddx,ddy
60  real :: rhoprof(2),rhoi
61  real :: xl,yl,wx,wy,w
62  real,parameter :: factor=.596831, hxmax=6.0, hymax=4.0, hzmax=150.
63
64
65  ! For forward simulations, make a loop over the number of species;
66  ! for backward simulations, make an additional loop over the
67  ! releasepoints
68  !***************************************************************************
69
70
71  do i=1,numpart
72!   if (itra1(i).ne.itime .or. abs(npoint(i)).gt.numpoint) goto 20
73    if (itra1(i).ne.itime ) goto 20
74
75  ! Determine age class of the particle
76    itage=abs(itra1(i)-itramem(i))
77    do nage=1,nageclass
78      if (itage.lt.lage(nage)) goto 33
79    end do
8033   continue
81
82
83  ! For special runs, interpolate the air density to the particle position
84  !************************************************************************
85  !***********************************************************************
86  !AF IND_SOURCE switches between different units for concentrations at the source
87  !Af    NOTE that in backward simulations the release of particles takes place
88  !Af    at the receptor and the sampling at the source.
89  !Af          1="mass"
90  !Af          2="mass mixing ratio"
91  !Af IND_RECEPTOR switches between different units for concentrations at the receptor
92  !Af          1="mass"
93  !Af          2="mass mixing ratio"
94
95  !Af switches for the conccalcfile:
96  !AF IND_SAMP =  0 : xmass * 1
97  !Af IND_SAMP = -1 : xmass / rho
98
99  !Af ind_samp is defined in readcommand.f
100
101    if ( ind_samp .eq. -1 ) then
102
103      ix=int(xtra1(i))
104      jy=int(ytra1(i))
105      ixp=ix+1
106      jyp=jy+1
107      ddx=xtra1(i)-real(ix)
108      ddy=ytra1(i)-real(jy)
109      rddx=1.-ddx
110      rddy=1.-ddy
111      p1=rddx*rddy
112      p2=ddx*rddy
113      p3=rddx*ddy
114      p4=ddx*ddy
115
116      do il=2,nz
117        if (height(il).gt.ztra1(i)) then
118          indz=il-1
119          indzp=il
120          goto 6
121        endif
122      end do
1236     continue
124
125      dz1=ztra1(i)-height(indz)
126      dz2=height(indzp)-ztra1(i)
127      dz=1./(dz1+dz2)
128
129  ! Take density from 2nd wind field in memory (accurate enough, no time interpolation needed)
130  !*****************************************************************************
131       if (xtra1(i).ne.xtra1(i))  print*,i,itra1(i),xtra1(i),ytra1(i),ix,jy
132      do ind=indz,indzp
133        rhoprof(ind-indz+1)=p1*rho(ix ,jy ,ind,2) &
134             +p2*rho(ixp,jy ,ind,2) &
135             +p3*rho(ix ,jyp,ind,2) &
136             +p4*rho(ixp,jyp,ind,2)
137      end do
138      rhoi=(dz1*rhoprof(2)+dz2*rhoprof(1))*dz
139   elseif (ind_samp.eq.0) then
140      rhoi = 1.
141   endif
142
143
144  !****************************************************************************
145  ! 1. Evaluate grid concentrations using a uniform kernel of bandwidths dx, dy
146  !****************************************************************************
147
148
149  ! For backward simulations, look from which release point the particle comes from
150  ! For domain-filling trajectory option, npoint contains a consecutive particle
151  ! number, not the release point information. Therefore, nrelpointer is set to 1
152  ! for the domain-filling option.
153  !*****************************************************************************
154
155    if ((ioutputforeachrelease.eq.0).or.(mdomainfill.ge.1)) then
156       nrelpointer=1
157    else
158       nrelpointer=npoint(i)
159    endif
160
161    do kz=1,numzgrid                ! determine height of cell
162      if (outheight(kz).gt.ztra1(i)) goto 21
163    end do
16421   continue
165    if (kz.le.numzgrid) then           ! inside output domain
166
167
168  !********************************
169  ! Do everything for mother domain
170  !********************************
171
172      xl=(xtra1(i)*dx+xoutshift)/dxout
173      yl=(ytra1(i)*dy+youtshift)/dyout
174      ix=int(xl)
175      if (xl.lt.0.) ix=ix-1
176      jy=int(yl)
177      if (yl.lt.0.) jy=jy-1
178
179  ! if (i.eq.10000) write(*,*) itime,xtra1(i),ytra1(i),ztra1(i),xl,yl
180
181
182  ! For particles aged less than 3 hours, attribute particle mass to grid cell
183  ! it resides in rather than use the kernel, in order to avoid its smoothing effect.
184  ! For older particles, use the uniform kernel.
185  ! If a particle is close to the domain boundary, do not use the kernel either.
186  !*****************************************************************************
187
188      if ((itage.lt.7200).or.(xl.lt.0.5).or.(yl.lt.0.5).or. &
189           (xl.gt.real(numxgrid-1)-0.5).or. &
190           (yl.gt.real(numygrid-1)-0.5)) then             ! no kernel, direct attribution to grid cell
191        if ((ix.ge.0).and.(jy.ge.0).and.(ix.le.numxgrid-1).and. &
192             (jy.le.numygrid-1)) then
193          do ks=1,nspec
194!        print*,'bef ks'
195!         print*,ix,jy,kz,ks,nrelpointer,nclass(i),nage
196          incl=nclass(i)
197             gridunc(ix,jy,kz,ks,nrelpointer,incl,nage)= &
198                  gridunc(ix,jy,kz,ks,nrelpointer,incl,nage)+ &
199                  xmass1(i,ks)/rhoi*weight
200          end do
201        endif
202
203      else                                 ! attribution via uniform kernel
204
205        ddx=xl-real(ix)                   ! distance to left cell border
206        ddy=yl-real(jy)                   ! distance to lower cell border
207        if (ddx.gt.0.5) then
208          ixp=ix+1
209          wx=1.5-ddx
210        else
211          ixp=ix-1
212          wx=0.5+ddx
213        endif
214
215        if (ddy.gt.0.5) then
216          jyp=jy+1
217          wy=1.5-ddy
218        else
219          jyp=jy-1
220          wy=0.5+ddy
221        endif
222
223
224  ! Determine mass fractions for four grid points
225  !**********************************************
226
227        if ((ix.ge.0).and.(ix.le.numxgrid-1)) then
228          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
229            w=wx*wy
230            do ks=1,nspec
231              gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
232                   gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
233                   xmass1(i,ks)/rhoi*weight*w
234            end do
235          endif
236
237          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
238            w=wx*(1.-wy)
239            do ks=1,nspec
240              gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
241                   gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
242                   xmass1(i,ks)/rhoi*weight*w
243            end do
244          endif
245        endif
246
247
248        if ((ixp.ge.0).and.(ixp.le.numxgrid-1)) then
249          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
250            w=(1.-wx)*(1.-wy)
251            do ks=1,nspec
252              gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
253                   gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
254                   xmass1(i,ks)/rhoi*weight*w
255            end do
256          endif
257
258          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
259            w=(1.-wx)*wy
260            do ks=1,nspec
261              gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
262                   gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
263                   xmass1(i,ks)/rhoi*weight*w
264            end do
265          endif
266        endif
267      endif
268
269
270
271  !************************************
272  ! Do everything for the nested domain
273  !************************************
274
275      if (nested_output.eq.1) then
276        xl=(xtra1(i)*dx+xoutshiftn)/dxoutn
277        yl=(ytra1(i)*dy+youtshiftn)/dyoutn
278        ix=int(xl)
279        if (xl.lt.0.) ix=ix-1
280        jy=int(yl)
281        if (yl.lt.0.) jy=jy-1
282
283
284  ! For particles aged less than 3 hours, attribute particle mass to grid cell
285  ! it resides in rather than use the kernel, in order to avoid its smoothing effect.
286  ! For older particles, use the uniform kernel.
287  ! If a particle is close to the domain boundary, do not use the kernel either.
288  !*****************************************************************************
289
290        if ((itage.lt.7200).or.(xl.lt.0.5).or.(yl.lt.0.5).or. &
291             (xl.gt.real(numxgridn-1)-0.5).or. &
292             (yl.gt.real(numygridn-1)-0.5)) then             ! no kernel, direct attribution to grid cell
293          if ((ix.ge.0).and.(jy.ge.0).and.(ix.le.numxgridn-1).and. &
294               (jy.le.numygridn-1)) then
295            do ks=1,nspec
296              griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
297                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
298                   xmass1(i,ks)/rhoi*weight
299            end do
300          endif
301
302        else                                 ! attribution via uniform kernel
303
304          ddx=xl-real(ix)                   ! distance to left cell border
305          ddy=yl-real(jy)                   ! distance to lower cell border
306          if (ddx.gt.0.5) then
307            ixp=ix+1
308            wx=1.5-ddx
309          else
310            ixp=ix-1
311            wx=0.5+ddx
312          endif
313
314          if (ddy.gt.0.5) then
315            jyp=jy+1
316            wy=1.5-ddy
317          else
318            jyp=jy-1
319            wy=0.5+ddy
320          endif
321
322
323  ! Determine mass fractions for four grid points
324  !**********************************************
325
326          if ((ix.ge.0).and.(ix.le.numxgridn-1)) then
327            if ((jy.ge.0).and.(jy.le.numygridn-1)) then
328              w=wx*wy
329              do ks=1,nspec
330                griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
331                     griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
332                     xmass1(i,ks)/rhoi*weight*w
333              end do
334            endif
335
336            if ((jyp.ge.0).and.(jyp.le.numygridn-1)) then
337              w=wx*(1.-wy)
338              do ks=1,nspec
339                griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
340                     griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
341                     xmass1(i,ks)/rhoi*weight*w
342              end do
343            endif
344          endif
345
346
347          if ((ixp.ge.0).and.(ixp.le.numxgridn-1)) then
348            if ((jyp.ge.0).and.(jyp.le.numygridn-1)) then
349              w=(1.-wx)*(1.-wy)
350              do ks=1,nspec
351                griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
352                     griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
353                     xmass1(i,ks)/rhoi*weight*w
354              end do
355            endif
356
357            if ((jy.ge.0).and.(jy.le.numygridn-1)) then
358              w=(1.-wx)*wy
359              do ks=1,nspec
360                griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
361                     griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
362                     xmass1(i,ks)/rhoi*weight*w
363              end do
364            endif
365          endif
366        endif
367
368      endif
369    endif
37020  continue
371  end do
372
373  !***********************************************************************
374  ! 2. Evaluate concentrations at receptor points, using the kernel method
375  !***********************************************************************
376
377  do n=1,numreceptor
378
379
380  ! Reset concentrations
381  !*********************
382
383    do ks=1,nspec
384      c(ks)=0.
385    end do
386
387
388  ! Estimate concentration at receptor
389  !***********************************
390
391    do i=1,numpart
392
393      if (itra1(i).ne.itime) goto 40
394      itage=abs(itra1(i)-itramem(i))
395
396      hz=min(50.+0.3*sqrt(real(itage)),hzmax)
397      zd=ztra1(i)/hz
398      if (zd.gt.1.) goto 40          ! save computing time, leave loop
399
400      hx=min((0.29+2.222e-3*sqrt(real(itage)))*dx+ &
401           real(itage)*1.2e-5,hxmax)                     ! 80 km/day
402      xd=(xtra1(i)-xreceptor(n))/hx
403      if (xd*xd.gt.1.) goto 40       ! save computing time, leave loop
404
405      hy=min((0.18+1.389e-3*sqrt(real(itage)))*dy+ &
406           real(itage)*7.5e-6,hymax)                     ! 80 km/day
407      yd=(ytra1(i)-yreceptor(n))/hy
408      if (yd*yd.gt.1.) goto 40       ! save computing time, leave loop
409      h=hx*hy*hz
410
411      r2=xd*xd+yd*yd+zd*zd
412      if (r2.lt.1.) then
413        xkern=factor*(1.-r2)
414        do ks=1,nspec
415          c(ks)=c(ks)+xmass1(i,ks)*xkern/h
416        end do
417      endif
41840    continue
419    end do
420
421    do ks=1,nspec
422      creceptor(n,ks)=creceptor(n,ks)+2.*weight*c(ks)/receptorarea(n)
423    end do
424  end do
425
426end subroutine conccalc_irreg
hosted by ZAMG