source: flexpart.git/src/conccalc.f90 @ 16b61a5

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

Reworked the domain-filling option (MPI). Fixed a slow loop which had errors in loop counter (MPI)

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