source: flexpart.git/src/conccalc_mpi.f90 @ 54cbd6c

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 54cbd6c was 5f9d14a, checked in by Espen Sollum ATMOS <eso@…>, 9 years ago

Updated wet depo scheme

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