source: flexpart.git/src/conccalc_mpi.f90 @ 47f96e5

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

Added (hardcoded) option to output number of particles per grid cell.

  • Property mode set to 100644
File size: 19.5 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          if (DRYBKDEP.or.WETBKDEP) then
208            do ks=1,nspec
209              gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
210                   gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
211                   xmass1(i,ks)/rhoi*weight*max(xscav_frac1(i,ks),0.0)
212            end do
213          else
214            if (lparticlecountoutput) then
215              do ks=1,nspec
216                gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
217                     gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+1
218              end do
219            else
220              do ks=1,nspec
221                gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
222                     gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
223                     xmass1(i,ks)/rhoi*weight
224              end do
225            end if
226          endif
227        endif
228
229      else                                 ! attribution via uniform kernel
230
231        ddx=xl-real(ix)                   ! distance to left cell border
232        ddy=yl-real(jy)                   ! distance to lower cell border
233        if (ddx.gt.0.5) then
234          ixp=ix+1
235          wx=1.5-ddx
236        else
237          ixp=ix-1
238          wx=0.5+ddx
239        endif
240
241        if (ddy.gt.0.5) then
242          jyp=jy+1
243          wy=1.5-ddy
244        else
245          jyp=jy-1
246          wy=0.5+ddy
247        endif
248
249
250  ! Determine mass fractions for four grid points
251  !**********************************************
252
253        if ((ix.ge.0).and.(ix.le.numxgrid-1)) then
254          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
255            w=wx*wy
256            if (DRYBKDEP.or.WETBKDEP) then
257               do ks=1,nspec
258                 gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
259                   gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
260                   xmass1(i,ks)/rhoi*w*weight*max(xscav_frac1(i,ks),0.0)
261               end do
262            else
263               do ks=1,nspec
264                 gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
265                   gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
266                   xmass1(i,ks)/rhoi*weight*w
267               end do
268            endif
269          endif
270
271          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
272            w=wx*(1.-wy)
273            if (DRYBKDEP.or.WETBKDEP) then
274              do ks=1,nspec
275                 gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
276                   gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
277                   xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
278               end do
279             else
280              do ks=1,nspec
281                 gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
282                   gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
283                   xmass1(i,ks)/rhoi*weight*w
284               end do
285             endif
286          endif
287        endif
288
289
290        if ((ixp.ge.0).and.(ixp.le.numxgrid-1)) then
291          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
292            w=(1.-wx)*(1.-wy)
293            if (DRYBKDEP.or.WETBKDEP) then
294               do ks=1,nspec
295                 gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
296                   gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
297                   xmass1(i,ks)/rhoi*w*weight*max(xscav_frac1(i,ks),0.0)
298               end do
299            else
300               do ks=1,nspec
301                 gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
302                   gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
303                   xmass1(i,ks)/rhoi*weight*w
304               end do
305            endif
306          endif
307
308          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
309            w=(1.-wx)*wy
310            if (DRYBKDEP.or.WETBKDEP) then
311               do ks=1,nspec
312                 gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
313                   gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
314                   xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
315               end do
316            else
317               do ks=1,nspec
318                 gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
319                   gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
320                   xmass1(i,ks)/rhoi*weight*w
321               end do
322            endif
323          endif
324        endif
325      endif
326
327
328
329  !************************************
330  ! Do everything for the nested domain
331  !************************************
332
333      if (nested_output.eq.1) then
334        xl=(xtra1(i)*dx+xoutshiftn)/dxoutn
335        yl=(ytra1(i)*dy+youtshiftn)/dyoutn
336        ix=int(xl)
337        if (xl.lt.0.) ix=ix-1
338        jy=int(yl)
339        if (yl.lt.0.) jy=jy-1
340
341
342  ! For particles aged less than 3 hours, attribute particle mass to grid cell
343  ! it resides in rather than use the kernel, in order to avoid its smoothing effect.
344  ! For older particles, use the uniform kernel.
345  ! If a particle is close to the domain boundary, do not use the kernel either.
346  !*****************************************************************************
347
348        if ((itage.lt.10800).or.(xl.lt.0.5).or.(yl.lt.0.5).or. &
349             (xl.gt.real(numxgridn-1)-0.5).or. &
350             (yl.gt.real(numygridn-1)-0.5)) then             ! no kernel, direct attribution to grid cell
351          if ((ix.ge.0).and.(jy.ge.0).and.(ix.le.numxgridn-1).and. &
352               (jy.le.numygridn-1)) then
353            if (DRYBKDEP.or.WETBKDEP) then
354               do ks=1,nspec
355                 griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
356                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
357                   xmass1(i,ks)/rhoi*weight*max(xscav_frac1(i,ks),0.0)
358               end do
359            else
360              if (lparticlecountoutput) then
361                do ks=1,nspec
362                  griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
363                       griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+1
364                end do
365              else
366                do ks=1,nspec
367                  griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
368                       griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
369                       xmass1(i,ks)/rhoi*weight
370                end do
371              endif
372            endif
373          endif
374
375        else                                 ! attribution via uniform kernel
376
377          ddx=xl-real(ix)                   ! distance to left cell border
378          ddy=yl-real(jy)                   ! distance to lower cell border
379          if (ddx.gt.0.5) then
380            ixp=ix+1
381            wx=1.5-ddx
382          else
383            ixp=ix-1
384            wx=0.5+ddx
385          endif
386
387          if (ddy.gt.0.5) then
388            jyp=jy+1
389            wy=1.5-ddy
390          else
391            jyp=jy-1
392            wy=0.5+ddy
393          endif
394
395
396  ! Determine mass fractions for four grid points
397  !**********************************************
398
399          if ((ix.ge.0).and.(ix.le.numxgridn-1)) then
400            if ((jy.ge.0).and.(jy.le.numygridn-1)) then
401              w=wx*wy
402              if (DRYBKDEP.or.WETBKDEP) then
403                 do ks=1,nspec
404                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
405                     griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
406                     xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
407                 end do
408              else
409                do ks=1,nspec
410                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
411                     griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
412                     xmass1(i,ks)/rhoi*weight*w
413                 end do
414              endif
415            endif
416
417            if ((jyp.ge.0).and.(jyp.le.numygridn-1)) then
418              w=wx*(1.-wy)
419              if (DRYBKDEP.or.WETBKDEP) then
420                 do ks=1,nspec
421                   griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
422                     griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
423                     xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
424                 end do
425              else
426                 do ks=1,nspec
427                   griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
428                     griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
429                     xmass1(i,ks)/rhoi*weight*w
430                 end do
431              endif
432            endif
433          endif
434
435
436          if ((ixp.ge.0).and.(ixp.le.numxgridn-1)) then
437            if ((jyp.ge.0).and.(jyp.le.numygridn-1)) then
438              w=(1.-wx)*(1.-wy)
439              if (DRYBKDEP.or.WETBKDEP) then
440                 do ks=1,nspec
441                   griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
442                     griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
443                     xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
444                 end do
445              else
446                 do ks=1,nspec
447                   griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
448                     griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
449                     xmass1(i,ks)/rhoi*weight*w
450                 end do
451              endif
452            endif
453
454            if ((jy.ge.0).and.(jy.le.numygridn-1)) then
455              w=(1.-wx)*wy
456              if (DRYBKDEP.or.WETBKDEP) then
457                 do ks=1,nspec
458                   griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
459                     griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
460                     xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
461                 end do
462              else
463                 do ks=1,nspec
464                    griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
465                     griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
466                     xmass1(i,ks)/rhoi*weight*w
467                 end do
468              endif
469            endif
470          endif
471        endif
472
473      endif
474    endif
47520  continue
476  end do
477
478  !***********************************************************************
479  ! 2. Evaluate concentrations at receptor points, using the kernel method
480  !***********************************************************************
481
482  do n=1,numreceptor
483
484
485  ! Reset concentrations
486  !*********************
487
488    do ks=1,nspec
489      c(ks)=0.
490    end do
491
492
493  ! Estimate concentration at receptor
494  !***********************************
495
496    do i=1,numpart
497
498      if (itra1(i).ne.itime) goto 40
499      itage=abs(itra1(i)-itramem(i))
500
501      hz=min(50.+0.3*sqrt(real(itage)),hzmax)
502      zd=ztra1(i)/hz
503      if (zd.gt.1.) goto 40          ! save computing time, leave loop
504
505      hx=min((0.29+2.222e-3*sqrt(real(itage)))*dx+ &
506           real(itage)*1.2e-5,hxmax)                     ! 80 km/day
507      xd=(xtra1(i)-xreceptor(n))/hx
508      if (xd*xd.gt.1.) goto 40       ! save computing time, leave loop
509
510      hy=min((0.18+1.389e-3*sqrt(real(itage)))*dy+ &
511           real(itage)*7.5e-6,hymax)                     ! 80 km/day
512      yd=(ytra1(i)-yreceptor(n))/hy
513      if (yd*yd.gt.1.) goto 40       ! save computing time, leave loop
514      h=hx*hy*hz
515
516      r2=xd*xd+yd*yd+zd*zd
517      if (r2.lt.1.) then
518        xkern=factor*(1.-r2)
519        do ks=1,nspec
520          c(ks)=c(ks)+xmass1(i,ks)*xkern/h
521        end do
522      endif
52340    continue
524    end do
525
526    do ks=1,nspec
527      creceptor(n,ks)=creceptor(n,ks)+2.*weight*c(ks)/receptorarea(n)
528    end do
529  end do
530
531  if (mp_measure_time) call mpif_mtime('conccalc',1)
532
533end subroutine conccalc
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG