source: flexpart.git/src/conccalc_mpi.f90 @ e52967c

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

Quick fix for segmentation fault when particle is at north pole

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