source: branches/flexpart91_hasod/src_parallel/conccalc.f90 @ 8

Last change on this file since 8 was 8, checked in by hasod, 11 years ago

Added parallel version of Flexpart91

File size: 15.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#if (defined _OPENMP)
52  use omp_lib
53#endif
54
55  implicit none
56
57  integer :: itime,itage,i,ix,jy,ixp,jyp,kz,ks,n,nage
58  integer :: il,ind,indz,indzp,nrelpointer
59  real :: rddx,rddy,p1,p2,p3,p4,dz1,dz2,dz
60  real :: weight,hx,hy,hz,h,xd,yd,zd,xkern,r2,c(maxspec),ddx,ddy
61  real :: rhoprof(2),rhoi
62  real :: xl,yl,wx,wy,w
63  integer :: griclunc
64#if (defined _OPENMP)
65  integer :: thread
66#endif
67  real,parameter :: factor=.596831, hxmax=6.0, hymax=4.0, hzmax=150.
68
69
70  ! For forward simulations, make a loop over the number of species;
71  ! for backward simulations, make an additional loop over the
72  ! releasepoints
73  !***************************************************************************
74
75!$OMP  PARALLEL PRIVATE(griclunc, itage,ix,jy,ixp,jyp,ddx,ddy,rddx,rddy, &
76!$OMP    p1,p2,p3,p4,nage,il,indz,indzp,dz1,dz2,dz,rhoprof,rhoi,      &
77!$OMP    thread,nrelpointer,kz,xl,yl,wx,wy,  &
78!$OMP    ks,i)
79
80#if (defined _OPENMP)
81      thread = OMP_GET_THREAD_NUM()
82#endif
83
84!$OMP DO
85
86  particleloop: do i=1,numpart
87
88    if (itra1(i).ne.itime) cycle
89
90    griclunc = mod(nclass(i)-1, nclassunc) + 1
91
92
93  ! Determine age class of the particle
94    itage=abs(itra1(i)-itramem(i))
95    do nage=1,nageclass
96      if (itage.lt.lage(nage)) goto 33
97    end do
9833   continue
99
100
101  ! For special runs, interpolate the air density to the particle position
102  !************************************************************************
103  !***********************************************************************
104  !AF IND_SOURCE switches between different units for concentrations at the source
105  !Af    NOTE that in backward simulations the release of particles takes place
106  !Af    at the receptor and the sampling at the source.
107  !Af          1="mass"
108  !Af          2="mass mixing ratio"
109  !Af IND_RECEPTOR switches between different units for concentrations at the receptor
110  !Af          1="mass"
111  !Af          2="mass mixing ratio"
112
113  !Af switches for the conccalcfile:
114  !AF IND_SAMP =  0 : xmass * 1
115  !Af IND_SAMP = -1 : xmass / rho
116
117  !Af ind_samp is defined in readcommand.f
118
119    if ( ind_samp .eq. -1 ) then
120
121      ix=int(xtra1(i))
122      jy=int(ytra1(i))
123      ixp=ix+1
124      jyp=jy+1
125      ddx=xtra1(i)-real(ix)
126      ddy=ytra1(i)-real(jy)
127      rddx=1.-ddx
128      rddy=1.-ddy
129      p1=rddx*rddy
130      p2=ddx*rddy
131      p3=rddx*ddy
132      p4=ddx*ddy
133
134      do il=2,nz
135        if (height(il).gt.ztra1(i)) then
136          indz=il-1
137          indzp=il
138          goto 6
139        endif
140      end do
1416     continue
142
143      dz1=ztra1(i)-height(indz)
144      dz2=height(indzp)-ztra1(i)
145      dz=1./(dz1+dz2)
146
147  ! Take density from 2nd wind field in memory (accurate enough, no time interpolation needed)
148  !*****************************************************************************
149      do ind=indz,indzp
150        rhoprof(ind-indz+1)=p1*rho(ix ,jy ,ind,2) &
151             +p2*rho(ixp,jy ,ind,2) &
152             +p3*rho(ix ,jyp,ind,2) &
153             +p4*rho(ixp,jyp,ind,2)
154      end do
155      rhoi=(dz1*rhoprof(2)+dz2*rhoprof(1))*dz
156   elseif (ind_samp.eq.0) then
157      rhoi = 1.
158   endif
159
160
161  !****************************************************************************
162  ! 1. Evaluate grid concentrations using a uniform kernel of bandwidths dx, dy
163  !****************************************************************************
164
165
166  ! For backward simulations, look from which release point the particle comes from
167  ! For domain-filling trajectory option, npoint contains a consecutive particle
168  ! number, not the release point information. Therefore, nrelpointer is set to 1
169  ! for the domain-filling option.
170  !*****************************************************************************
171
172    if ((ioutputforeachrelease.eq.0).or.(mdomainfill.eq.1)) then
173       nrelpointer=1
174    else
175       nrelpointer=npoint(i)
176    endif
177
178    do kz=1,numzgrid                ! determine height of cell
179      if (outheight(kz).gt.ztra1(i)) exit
180    end do
181
182    if (kz.le.numzgrid) then           ! inside output domain
183
184
185  !********************************
186  ! Do everything for mother domain
187  !********************************
188
189      xl=(xtra1(i)*dx+xoutshift)/dxout
190      yl=(ytra1(i)*dy+youtshift)/dyout
191      ix=int(xl)
192      if (xl.lt.0.) ix=ix-1
193      jy=int(yl)
194      if (yl.lt.0.) jy=jy-1
195
196  ! if (i.eq.10000) write(*,*) itime,xtra1(i),ytra1(i),ztra1(i),xl,yl
197
198
199  ! For particles aged less than 3 hours, attribute particle mass to grid cell
200  ! it resides in rather than use the kernel, in order to avoid its smoothing effect.
201  ! For older particles, use the uniform kernel.
202  ! If a particle is close to the domain boundary, do not use the kernel either.
203  !*****************************************************************************
204
205      if ((itage.lt.10800).or.(xl.lt.0.5).or.(yl.lt.0.5).or. &
206           (xl.gt.real(numxgrid-1)-0.5).or. &
207           (yl.gt.real(numygrid-1)-0.5)) then             ! no kernel, direct attribution to grid cell
208        if ((ix.ge.0).and.(jy.ge.0).and.(ix.le.numxgrid-1).and. &
209             (jy.le.numygrid-1)) then
210!$OMP CRITICAL
211          do ks=1,nspec
212            gridunc(ix,jy,kz,ks,nrelpointer,griclunc,nage)= &
213                 gridunc(ix,jy,kz,ks,nrelpointer,griclunc,nage)+ &
214                 xmass1(i,ks)/rhoi*weight
215          end do
216!$OMP END CRITICAL
217        endif
218
219      else                                 ! attribution via uniform kernel
220
221        ddx=xl-real(ix)                   ! distance to left cell border
222        ddy=yl-real(jy)                   ! distance to lower cell border
223        if (ddx.gt.0.5) then
224          ixp=ix+1
225          wx=1.5-ddx
226        else
227          ixp=ix-1
228          wx=0.5+ddx
229        endif
230
231        if (ddy.gt.0.5) then
232          jyp=jy+1
233          wy=1.5-ddy
234        else
235          jyp=jy-1
236          wy=0.5+ddy
237        endif
238
239
240  ! Determine mass fractions for four grid points
241  !**********************************************
242
243        if ((ix.ge.0).and.(ix.le.numxgrid-1)) then
244          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
245            w=wx*wy
246            do ks=1,nspec
247              gridunc(ix,jy,kz,ks,nrelpointer,griclunc,nage)= &
248                   gridunc(ix,jy,kz,ks,nrelpointer,griclunc,nage)+ &
249                   xmass1(i,ks)/rhoi*weight*w
250            end do
251          endif
252
253          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
254            w=wx*(1.-wy)
255            do ks=1,nspec
256              gridunc(ix,jyp,kz,ks,nrelpointer,griclunc,nage)= &
257                   gridunc(ix,jyp,kz,ks,nrelpointer,griclunc,nage)+ &
258                   xmass1(i,ks)/rhoi*weight*w
259            end do
260          endif
261        endif
262
263
264        if ((ixp.ge.0).and.(ixp.le.numxgrid-1)) then
265          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
266            w=(1.-wx)*(1.-wy)
267            do ks=1,nspec
268              gridunc(ixp,jyp,kz,ks,nrelpointer,griclunc,nage)= &
269                   gridunc(ixp,jyp,kz,ks,nrelpointer,griclunc,nage)+ &
270                   xmass1(i,ks)/rhoi*weight*w
271            end do
272          endif
273
274          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
275            w=(1.-wx)*wy
276            do ks=1,nspec
277              gridunc(ixp,jy,kz,ks,nrelpointer,griclunc,nage)= &
278                   gridunc(ixp,jy,kz,ks,nrelpointer,griclunc,nage)+ &
279                   xmass1(i,ks)/rhoi*weight*w
280            end do
281          endif
282        endif
283      endif
284
285
286
287  !************************************
288  ! Do everything for the nested domain
289  !************************************
290
291      if (nested_output.eq.1) then
292        xl=(xtra1(i)*dx+xoutshiftn)/dxoutn
293        yl=(ytra1(i)*dy+youtshiftn)/dyoutn
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.10800).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,griclunc,nage)= &
313                   griduncn(ix,jy,kz,ks,nrelpointer,griclunc,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,griclunc,nage)= &
347                     griduncn(ix,jy,kz,ks,nrelpointer,griclunc,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,griclunc,nage)= &
356                     griduncn(ix,jyp,kz,ks,nrelpointer,griclunc,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,griclunc,nage)= &
368                     griduncn(ixp,jyp,kz,ks,nrelpointer,griclunc,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,griclunc,nage)= &
377                     griduncn(ixp,jy,kz,ks,nrelpointer,griclunc,nage)+ &
378                     xmass1(i,ks)/rhoi*weight*w
379              end do
380            endif
381          endif
382        endif
383
384      endif
385    endif
386
387  end do particleloop
388
389!$OMP END DO
390
391!$OMP END PARALLEL
392
393
394  !***********************************************************************
395  ! 2. Evaluate concentrations at receptor points, using the kernel method
396  !***********************************************************************
397
398  do n=1,numreceptor
399
400
401  ! Reset concentrations
402  !*********************
403
404    do ks=1,nspec
405      c(ks)=0.
406    end do
407
408
409  ! Estimate concentration at receptor
410  !***********************************
411
412    do i=1,numpart
413
414      if (itra1(i).ne.itime) goto 40
415      itage=abs(itra1(i)-itramem(i))
416
417      hz=min(50.+0.3*sqrt(real(itage)),hzmax)
418      zd=ztra1(i)/hz
419      if (zd.gt.1.) goto 40          ! save computing time, leave loop
420
421      hx=min((0.29+2.222e-3*sqrt(real(itage)))*dx+ &
422           real(itage)*1.2e-5,hxmax)                     ! 80 km/day
423      xd=(xtra1(i)-xreceptor(n))/hx
424      if (xd*xd.gt.1.) goto 40       ! save computing time, leave loop
425
426      hy=min((0.18+1.389e-3*sqrt(real(itage)))*dy+ &
427           real(itage)*7.5e-6,hymax)                     ! 80 km/day
428      yd=(ytra1(i)-yreceptor(n))/hy
429      if (yd*yd.gt.1.) goto 40       ! save computing time, leave loop
430      h=hx*hy*hz
431
432      r2=xd*xd+yd*yd+zd*zd
433      if (r2.lt.1.) then
434        xkern=factor*(1.-r2)
435        do ks=1,nspec
436          c(ks)=c(ks)+xmass1(i,ks)*xkern/h
437        end do
438      endif
43940    continue
440    end do
441
442    do ks=1,nspec
443      creceptor(n,ks)=creceptor(n,ks)+2.*weight*c(ks)/receptorarea(n)
444    end do
445  end do
446
447end subroutine conccalc
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG