source: flexpart.git/src/conccalc_mpi.f90 @ 8a65cb0

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

Added code, makefile for dev branch

  • Property mode set to 100755
File size: 15.0 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
65
66  ! For forward simulations, make a loop over the number of species;
67  ! for backward simulations, make an additional loop over the
68  ! releasepoints
69  !***************************************************************************
70
71  if (mp_measure_time) call mpif_mtime('conccalc',0)
72
73  do i=1,numpart
74    if (itra1(i).ne.itime) goto 20
75
76  ! Determine age class of the particle
77    itage=abs(itra1(i)-itramem(i))
78    do nage=1,nageclass
79      if (itage.lt.lage(nage)) goto 33
80    end do
8133   continue
82
83
84  ! For special runs, interpolate the air density to the particle position
85  !************************************************************************
86  !***********************************************************************
87  !AF IND_SOURCE switches between different units for concentrations at the source
88  !Af    NOTE that in backward simulations the release of particles takes place
89  !Af    at the receptor and the sampling at the source.
90  !Af          1="mass"
91  !Af          2="mass mixing ratio"
92  !Af IND_RECEPTOR switches between different units for concentrations at the receptor
93  !Af          1="mass"
94  !Af          2="mass mixing ratio"
95
96  !Af switches for the conccalcfile:
97  !AF IND_SAMP =  0 : xmass * 1
98  !Af IND_SAMP = -1 : xmass / rho
99
100  !Af ind_samp is defined in readcommand.f
101
102    if ( ind_samp .eq. -1 ) then
103
104      ix=int(xtra1(i))
105      jy=int(ytra1(i))
106      ixp=ix+1
107      jyp=jy+1
108      ddx=xtra1(i)-real(ix)
109      ddy=ytra1(i)-real(jy)
110      rddx=1.-ddx
111      rddy=1.-ddy
112      p1=rddx*rddy
113      p2=ddx*rddy
114      p3=rddx*ddy
115      p4=ddx*ddy
116
117      do il=2,nz
118        if (height(il).gt.ztra1(i)) then
119          indz=il-1
120          indzp=il
121          goto 6
122        endif
123      end do
1246     continue
125
126      dz1=ztra1(i)-height(indz)
127      dz2=height(indzp)-ztra1(i)
128      dz=1./(dz1+dz2)
129
130  ! Take density from 2nd wind field in memory (accurate enough, no time interpolation needed)
131  !*****************************************************************************
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.eq.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.10800).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            gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
195                 gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
196                 xmass1(i,ks)/rhoi*weight
197          end do
198        endif
199
200      else                                 ! attribution via uniform kernel
201
202        ddx=xl-real(ix)                   ! distance to left cell border
203        ddy=yl-real(jy)                   ! distance to lower cell border
204        if (ddx.gt.0.5) then
205          ixp=ix+1
206          wx=1.5-ddx
207        else
208          ixp=ix-1
209          wx=0.5+ddx
210        endif
211
212        if (ddy.gt.0.5) then
213          jyp=jy+1
214          wy=1.5-ddy
215        else
216          jyp=jy-1
217          wy=0.5+ddy
218        endif
219
220
221  ! Determine mass fractions for four grid points
222  !**********************************************
223
224        if ((ix.ge.0).and.(ix.le.numxgrid-1)) then
225          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
226            w=wx*wy
227            do ks=1,nspec
228              gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
229                   gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
230                   xmass1(i,ks)/rhoi*weight*w
231            end do
232          endif
233
234          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
235            w=wx*(1.-wy)
236            do ks=1,nspec
237              gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
238                   gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
239                   xmass1(i,ks)/rhoi*weight*w
240            end do
241          endif
242        endif
243
244
245        if ((ixp.ge.0).and.(ixp.le.numxgrid-1)) then
246          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
247            w=(1.-wx)*(1.-wy)
248            do ks=1,nspec
249              gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
250                   gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
251                   xmass1(i,ks)/rhoi*weight*w
252            end do
253          endif
254
255          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
256            w=(1.-wx)*wy
257            do ks=1,nspec
258              gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
259                   gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
260                   xmass1(i,ks)/rhoi*weight*w
261            end do
262          endif
263        endif
264      endif
265
266
267
268  !************************************
269  ! Do everything for the nested domain
270  !************************************
271
272      if (nested_output.eq.1) then
273        xl=(xtra1(i)*dx+xoutshiftn)/dxoutn
274        yl=(ytra1(i)*dy+youtshiftn)/dyoutn
275        ix=int(xl)
276        if (xl.lt.0.) ix=ix-1
277        jy=int(yl)
278        if (yl.lt.0.) jy=jy-1
279
280
281  ! For particles aged less than 3 hours, attribute particle mass to grid cell
282  ! it resides in rather than use the kernel, in order to avoid its smoothing effect.
283  ! For older particles, use the uniform kernel.
284  ! If a particle is close to the domain boundary, do not use the kernel either.
285  !*****************************************************************************
286
287        if ((itage.lt.10800).or.(xl.lt.0.5).or.(yl.lt.0.5).or. &
288             (xl.gt.real(numxgridn-1)-0.5).or. &
289             (yl.gt.real(numygridn-1)-0.5)) then             ! no kernel, direct attribution to grid cell
290          if ((ix.ge.0).and.(jy.ge.0).and.(ix.le.numxgridn-1).and. &
291               (jy.le.numygridn-1)) then
292            do ks=1,nspec
293              griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
294                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
295                   xmass1(i,ks)/rhoi*weight
296            end do
297          endif
298
299        else                                 ! attribution via uniform kernel
300
301          ddx=xl-real(ix)                   ! distance to left cell border
302          ddy=yl-real(jy)                   ! distance to lower cell border
303          if (ddx.gt.0.5) then
304            ixp=ix+1
305            wx=1.5-ddx
306          else
307            ixp=ix-1
308            wx=0.5+ddx
309          endif
310
311          if (ddy.gt.0.5) then
312            jyp=jy+1
313            wy=1.5-ddy
314          else
315            jyp=jy-1
316            wy=0.5+ddy
317          endif
318
319
320  ! Determine mass fractions for four grid points
321  !**********************************************
322
323          if ((ix.ge.0).and.(ix.le.numxgridn-1)) then
324            if ((jy.ge.0).and.(jy.le.numygridn-1)) then
325              w=wx*wy
326              do ks=1,nspec
327                griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
328                     griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
329                     xmass1(i,ks)/rhoi*weight*w
330              end do
331            endif
332
333            if ((jyp.ge.0).and.(jyp.le.numygridn-1)) then
334              w=wx*(1.-wy)
335              do ks=1,nspec
336                griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
337                     griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
338                     xmass1(i,ks)/rhoi*weight*w
339              end do
340            endif
341          endif
342
343
344          if ((ixp.ge.0).and.(ixp.le.numxgridn-1)) then
345            if ((jyp.ge.0).and.(jyp.le.numygridn-1)) then
346              w=(1.-wx)*(1.-wy)
347              do ks=1,nspec
348                griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
349                     griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
350                     xmass1(i,ks)/rhoi*weight*w
351              end do
352            endif
353
354            if ((jy.ge.0).and.(jy.le.numygridn-1)) then
355              w=(1.-wx)*wy
356              do ks=1,nspec
357                griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
358                     griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
359                     xmass1(i,ks)/rhoi*weight*w
360              end do
361            endif
362          endif
363        endif
364
365      endif
366    endif
36720  continue
368  end do
369
370  !***********************************************************************
371  ! 2. Evaluate concentrations at receptor points, using the kernel method
372  !***********************************************************************
373
374  do n=1,numreceptor
375
376
377  ! Reset concentrations
378  !*********************
379
380    do ks=1,nspec
381      c(ks)=0.
382    end do
383
384
385  ! Estimate concentration at receptor
386  !***********************************
387
388    do i=1,numpart
389
390      if (itra1(i).ne.itime) goto 40
391      itage=abs(itra1(i)-itramem(i))
392
393      hz=min(50.+0.3*sqrt(real(itage)),hzmax)
394      zd=ztra1(i)/hz
395      if (zd.gt.1.) goto 40          ! save computing time, leave loop
396
397      hx=min((0.29+2.222e-3*sqrt(real(itage)))*dx+ &
398           real(itage)*1.2e-5,hxmax)                     ! 80 km/day
399      xd=(xtra1(i)-xreceptor(n))/hx
400      if (xd*xd.gt.1.) goto 40       ! save computing time, leave loop
401
402      hy=min((0.18+1.389e-3*sqrt(real(itage)))*dy+ &
403           real(itage)*7.5e-6,hymax)                     ! 80 km/day
404      yd=(ytra1(i)-yreceptor(n))/hy
405      if (yd*yd.gt.1.) goto 40       ! save computing time, leave loop
406      h=hx*hy*hz
407
408      r2=xd*xd+yd*yd+zd*zd
409      if (r2.lt.1.) then
410        xkern=factor*(1.-r2)
411        do ks=1,nspec
412          c(ks)=c(ks)+xmass1(i,ks)*xkern/h
413        end do
414      endif
41540    continue
416    end do
417
418    do ks=1,nspec
419      creceptor(n,ks)=creceptor(n,ks)+2.*weight*c(ks)/receptorarea(n)
420    end do
421  end do
422
423  if (mp_measure_time) call mpif_mtime('conccalc',1)
424
425end subroutine conccalc
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG