Ticket #74: conccalc_reg.2.f90

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