source: flexpart.git/src/conccalc_mpi.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 18.1 KB
Line 
1subroutine conccalc(itime,weight)
2  !                 i     i
3  !*****************************************************************************
4  !                                                                            *
5  !     Calculation of the concentrations on a regular grid using volume       *
6  !     sampling                                                               *
7  !                                                                            *
8  !     Author: A. Stohl                                                       *
9  !                                                                            *
10  !     24 May 1996                                                            *
11  !                                                                            *
12  !     April 2000: Update to calculate age spectra                            *
13  !                 Bug fix to avoid negative conc. at the domain boundaries,  *
14  !                 as suggested by Petra Seibert                              *
15  !                                                                            *
16  !     2 July 2002: re-order if-statements in order to optimize CPU time      *
17  !                                                                            *
18  !                                                                            *
19  !*****************************************************************************
20  !                                                                            *
21  ! Variables:                                                                 *
22  ! nspeciesdim     = nspec for forward runs, 1 for backward runs              *
23  !                                                                            *
24  !*****************************************************************************
25
26  use unc_mod
27  use outg_mod
28  use par_mod
29  use com_mod
30  use mpi_mod
31
32  implicit none
33
34  integer,intent(in) :: itime
35  integer :: itage,i,ix,jy,ixp,jyp,kz,ks,n,nage
36  integer :: il,ind,indz,indzp,nrelpointer
37  real,intent(in) :: weight
38  real :: rddx,rddy,p1,p2,p3,p4,dz1,dz2,dz
39  real :: hx,hy,hz,h,xd,yd,zd,xkern,r2,c(maxspec),ddx,ddy
40  real :: rhoprof(2),rhoi
41  real :: xl,yl,wx,wy,w
42  real,parameter :: factor=.596831, hxmax=6.0, hymax=4.0, hzmax=150.
43  integer :: mind2
44  ! mind2        eso: pointer to 2nd windfield in memory
45
46  ! For forward simulations, make a loop over the number of species;
47  ! for backward simulations, make an additional loop over the
48  ! releasepoints
49  !***************************************************************************
50
51  if (mp_measure_time) call mpif_mtime('conccalc',0)
52
53  mind2=memind(2)
54
55  do i=1,numpart
56    if (itra1(i).ne.itime) goto 20
57
58  ! Determine age class of the particle
59    itage=abs(itra1(i)-itramem(i))
60    do nage=1,nageclass
61      if (itage.lt.lage(nage)) goto 33
62    end do
6333   continue
64
65
66  ! For special runs, interpolate the air density to the particle position
67  !************************************************************************
68  !***********************************************************************
69  !AF IND_SOURCE switches between different units for concentrations at the source
70  !Af    NOTE that in backward simulations the release of particles takes place
71  !Af    at the receptor and the sampling at the source.
72  !Af          1="mass"
73  !Af          2="mass mixing ratio"
74  !Af IND_RECEPTOR switches between different units for concentrations at the receptor
75  !Af          1="mass"
76  !Af          2="mass mixing ratio"
77
78  !Af switches for the conccalcfile:
79  !AF IND_SAMP =  0 : xmass * 1
80  !Af IND_SAMP = -1 : xmass / rho
81
82  !Af ind_samp is defined in readcommand.f
83
84    if ( ind_samp .eq. -1 ) then
85
86      ix=int(xtra1(i))
87      jy=int(ytra1(i))
88      ixp=ix+1
89      jyp=jy+1
90      ddx=xtra1(i)-real(ix)
91      ddy=ytra1(i)-real(jy)
92      rddx=1.-ddx
93      rddy=1.-ddy
94      p1=rddx*rddy
95      p2=ddx*rddy
96      p3=rddx*ddy
97      p4=ddx*ddy
98
99! eso: Temporary fix for particle exactly at north pole
100      if (jyp >= nymax) then
101        ! write(*,*) 'WARNING: conccalc.f90 jyp >= nymax'
102        jyp=jyp-1
103      end if
104
105      do il=2,nz
106        if (height(il).gt.ztra1(i)) then
107          indz=il-1
108          indzp=il
109          goto 6
110        endif
111      end do
1126     continue
113
114      dz1=ztra1(i)-height(indz)
115      dz2=height(indzp)-ztra1(i)
116      dz=1./(dz1+dz2)
117
118  ! Take density from 2nd wind field in memory (accurate enough, no time interpolation needed)
119  !*****************************************************************************
120      do ind=indz,indzp
121        ! rhoprof(ind-indz+1)=p1*rho(ix ,jy ,ind,2) &
122        !      +p2*rho(ixp,jy ,ind,2) &
123        !      +p3*rho(ix ,jyp,ind,2) &
124        !      +p4*rho(ixp,jyp,ind,2)
125        rhoprof(ind-indz+1)=p1*rho(ix ,jy ,ind,mind2) &
126             +p2*rho(ixp,jy ,ind,mind2) &
127             +p3*rho(ix ,jyp,ind,mind2) &
128             +p4*rho(ixp,jyp,ind,mind2)
129
130      end do
131      rhoi=(dz1*rhoprof(2)+dz2*rhoprof(1))*dz
132   elseif (ind_samp.eq.0) then
133      rhoi = 1.
134   endif
135
136
137  !****************************************************************************
138  ! 1. Evaluate grid concentrations using a uniform kernel of bandwidths dx, dy
139  !****************************************************************************
140
141
142  ! For backward simulations, look from which release point the particle comes from
143  ! For domain-filling trajectory option, npoint contains a consecutive particle
144  ! number, not the release point information. Therefore, nrelpointer is set to 1
145  ! for the domain-filling option.
146  !*****************************************************************************
147
148    if ((ioutputforeachrelease.eq.0).or.(mdomainfill.eq.1)) then
149       nrelpointer=1
150    else
151       nrelpointer=npoint(i)
152    endif
153
154    do kz=1,numzgrid                ! determine height of cell
155      if (outheight(kz).gt.ztra1(i)) goto 21
156    end do
15721   continue
158    if (kz.le.numzgrid) then           ! inside output domain
159
160
161  !********************************
162  ! Do everything for mother domain
163  !********************************
164
165      xl=(xtra1(i)*dx+xoutshift)/dxout
166      yl=(ytra1(i)*dy+youtshift)/dyout
167      ix=int(xl)
168      if (xl.lt.0.) ix=ix-1
169      jy=int(yl)
170      if (yl.lt.0.) jy=jy-1
171
172  ! if (i.eq.10000) write(*,*) itime,xtra1(i),ytra1(i),ztra1(i),xl,yl
173
174
175  ! For particles aged less than 3 hours, attribute particle mass to grid cell
176  ! it resides in rather than use the kernel, in order to avoid its smoothing effect.
177  ! For older particles, use the uniform kernel.
178  ! If a particle is close to the domain boundary, do not use the kernel either.
179  !*****************************************************************************
180
181      if ((.not.lusekerneloutput).or.(itage.lt.10800).or. &
182           (xl.lt.0.5).or.(yl.lt.0.5).or. &
183           (xl.gt.real(numxgrid-1)-0.5).or. &
184           (yl.gt.real(numygrid-1)-0.5)) then             ! no kernel, direct attribution to grid cell
185        if ((ix.ge.0).and.(jy.ge.0).and.(ix.le.numxgrid-1).and. &
186             (jy.le.numygrid-1)) then
187          if (DRYBKDEP.or.WETBKDEP) then
188            do ks=1,nspec
189              gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
190                   gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
191                   xmass1(i,ks)/rhoi*weight*max(xscav_frac1(i,ks),0.0)
192            end do
193          else
194            if (lparticlecountoutput) then
195              do ks=1,nspec
196                gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
197                     gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+1
198              end do
199            else
200              do ks=1,nspec
201                gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
202                     gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
203                     xmass1(i,ks)/rhoi*weight
204              end do
205            end if
206          endif
207        endif
208
209      else                                 ! attribution via uniform kernel
210
211        ddx=xl-real(ix)                   ! distance to left cell border
212        ddy=yl-real(jy)                   ! distance to lower cell border
213        if (ddx.gt.0.5) then
214          ixp=ix+1
215          wx=1.5-ddx
216        else
217          ixp=ix-1
218          wx=0.5+ddx
219        endif
220
221        if (ddy.gt.0.5) then
222          jyp=jy+1
223          wy=1.5-ddy
224        else
225          jyp=jy-1
226          wy=0.5+ddy
227        endif
228
229
230  ! Determine mass fractions for four grid points
231  !**********************************************
232
233        if ((ix.ge.0).and.(ix.le.numxgrid-1)) then
234          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
235            w=wx*wy
236            if (DRYBKDEP.or.WETBKDEP) then
237               do ks=1,nspec
238                 gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
239                   gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
240                   xmass1(i,ks)/rhoi*w*weight*max(xscav_frac1(i,ks),0.0)
241               end do
242            else
243               do ks=1,nspec
244                 gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
245                   gridunc(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
246                   xmass1(i,ks)/rhoi*weight*w
247               end do
248            endif
249          endif
250
251          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
252            w=wx*(1.-wy)
253            if (DRYBKDEP.or.WETBKDEP) then
254              do ks=1,nspec
255                 gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
256                   gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
257                   xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
258               end do
259             else
260              do ks=1,nspec
261                 gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
262                   gridunc(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
263                   xmass1(i,ks)/rhoi*weight*w
264               end do
265             endif
266          endif
267        endif
268
269
270        if ((ixp.ge.0).and.(ixp.le.numxgrid-1)) then
271          if ((jyp.ge.0).and.(jyp.le.numygrid-1)) then
272            w=(1.-wx)*(1.-wy)
273            if (DRYBKDEP.or.WETBKDEP) then
274               do ks=1,nspec
275                 gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
276                   gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
277                   xmass1(i,ks)/rhoi*w*weight*max(xscav_frac1(i,ks),0.0)
278               end do
279            else
280               do ks=1,nspec
281                 gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
282                   gridunc(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
283                   xmass1(i,ks)/rhoi*weight*w
284               end do
285            endif
286          endif
287
288          if ((jy.ge.0).and.(jy.le.numygrid-1)) then
289            w=(1.-wx)*wy
290            if (DRYBKDEP.or.WETBKDEP) then
291               do ks=1,nspec
292                 gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
293                   gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
294                   xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
295               end do
296            else
297               do ks=1,nspec
298                 gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
299                   gridunc(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
300                   xmass1(i,ks)/rhoi*weight*w
301               end do
302            endif
303          endif
304        endif
305      endif
306
307
308
309  !************************************
310  ! Do everything for the nested domain
311  !************************************
312
313      if (nested_output.eq.1) then
314        xl=(xtra1(i)*dx+xoutshiftn)/dxoutn
315        yl=(ytra1(i)*dy+youtshiftn)/dyoutn
316        ix=int(xl)
317        if (xl.lt.0.) ix=ix-1
318        jy=int(yl)
319        if (yl.lt.0.) jy=jy-1
320
321
322  ! For particles aged less than 3 hours, attribute particle mass to grid cell
323  ! it resides in rather than use the kernel, in order to avoid its smoothing effect.
324  ! For older particles, use the uniform kernel.
325  ! If a particle is close to the domain boundary, do not use the kernel either.
326  !*****************************************************************************
327
328        if ((itage.lt.10800).or.(xl.lt.0.5).or.(yl.lt.0.5).or. &
329             (xl.gt.real(numxgridn-1)-0.5).or. &
330             (yl.gt.real(numygridn-1)-0.5).or.((.not.lusekerneloutput))) then
331! no kernel, direct attribution to grid cell
332          if ((ix.ge.0).and.(jy.ge.0).and.(ix.le.numxgridn-1).and. &
333               (jy.le.numygridn-1)) then
334            if (DRYBKDEP.or.WETBKDEP) then
335               do ks=1,nspec
336                 griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
337                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
338                   xmass1(i,ks)/rhoi*weight*max(xscav_frac1(i,ks),0.0)
339               end do
340            else
341              if (lparticlecountoutput) then
342                do ks=1,nspec
343                  griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
344                       griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+1
345                end do
346              else
347                do ks=1,nspec
348                  griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
349                       griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
350                       xmass1(i,ks)/rhoi*weight
351                end do
352              endif
353            endif
354          endif
355
356        else                                 ! attribution via uniform kernel
357
358          ddx=xl-real(ix)                   ! distance to left cell border
359          ddy=yl-real(jy)                   ! distance to lower cell border
360          if (ddx.gt.0.5) then
361            ixp=ix+1
362            wx=1.5-ddx
363          else
364            ixp=ix-1
365            wx=0.5+ddx
366          endif
367
368          if (ddy.gt.0.5) then
369            jyp=jy+1
370            wy=1.5-ddy
371          else
372            jyp=jy-1
373            wy=0.5+ddy
374          endif
375
376
377  ! Determine mass fractions for four grid points
378  !**********************************************
379
380          if ((ix.ge.0).and.(ix.le.numxgridn-1)) then
381            if ((jy.ge.0).and.(jy.le.numygridn-1)) then
382              w=wx*wy
383              if (DRYBKDEP.or.WETBKDEP) then
384                 do ks=1,nspec
385                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
386                     griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
387                     xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
388                 end do
389              else
390                do ks=1,nspec
391                   griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)= &
392                     griduncn(ix,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
393                     xmass1(i,ks)/rhoi*weight*w
394                 end do
395              endif
396            endif
397
398            if ((jyp.ge.0).and.(jyp.le.numygridn-1)) then
399              w=wx*(1.-wy)
400              if (DRYBKDEP.or.WETBKDEP) then
401                 do ks=1,nspec
402                   griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
403                     griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
404                     xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
405                 end do
406              else
407                 do ks=1,nspec
408                   griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
409                     griduncn(ix,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
410                     xmass1(i,ks)/rhoi*weight*w
411                 end do
412              endif
413            endif
414          endif
415
416
417          if ((ixp.ge.0).and.(ixp.le.numxgridn-1)) then
418            if ((jyp.ge.0).and.(jyp.le.numygridn-1)) then
419              w=(1.-wx)*(1.-wy)
420              if (DRYBKDEP.or.WETBKDEP) then
421                 do ks=1,nspec
422                   griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
423                     griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
424                     xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
425                 end do
426              else
427                 do ks=1,nspec
428                   griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)= &
429                     griduncn(ixp,jyp,kz,ks,nrelpointer,nclass(i),nage)+ &
430                     xmass1(i,ks)/rhoi*weight*w
431                 end do
432              endif
433            endif
434
435            if ((jy.ge.0).and.(jy.le.numygridn-1)) then
436              w=(1.-wx)*wy
437              if (DRYBKDEP.or.WETBKDEP) then
438                 do ks=1,nspec
439                   griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
440                     griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
441                     xmass1(i,ks)/rhoi*weight*w*max(xscav_frac1(i,ks),0.0)
442                 end do
443              else
444                 do ks=1,nspec
445                    griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)= &
446                     griduncn(ixp,jy,kz,ks,nrelpointer,nclass(i),nage)+ &
447                     xmass1(i,ks)/rhoi*weight*w
448                 end do
449              endif
450            endif
451          endif
452        endif
453
454      endif
455    endif
45620  continue
457  end do
458
459  !***********************************************************************
460  ! 2. Evaluate concentrations at receptor points, using the kernel method
461  !***********************************************************************
462
463  do n=1,numreceptor
464
465
466  ! Reset concentrations
467  !*********************
468
469    do ks=1,nspec
470      c(ks)=0.
471    end do
472
473
474  ! Estimate concentration at receptor
475  !***********************************
476
477    do i=1,numpart
478
479      if (itra1(i).ne.itime) goto 40
480      itage=abs(itra1(i)-itramem(i))
481
482      hz=min(50.+0.3*sqrt(real(itage)),hzmax)
483      zd=ztra1(i)/hz
484      if (zd.gt.1.) goto 40          ! save computing time, leave loop
485
486      hx=min((0.29+2.222e-3*sqrt(real(itage)))*dx+ &
487           real(itage)*1.2e-5,hxmax)                     ! 80 km/day
488      xd=(xtra1(i)-xreceptor(n))/hx
489      if (xd*xd.gt.1.) goto 40       ! save computing time, leave loop
490
491      hy=min((0.18+1.389e-3*sqrt(real(itage)))*dy+ &
492           real(itage)*7.5e-6,hymax)                     ! 80 km/day
493      yd=(ytra1(i)-yreceptor(n))/hy
494      if (yd*yd.gt.1.) goto 40       ! save computing time, leave loop
495      h=hx*hy*hz
496
497      r2=xd*xd+yd*yd+zd*zd
498      if (r2.lt.1.) then
499        xkern=factor*(1.-r2)
500        do ks=1,nspec
501          c(ks)=c(ks)+xmass1(i,ks)*xkern/h
502        end do
503      endif
50440    continue
505    end do
506
507    do ks=1,nspec
508      creceptor(n,ks)=creceptor(n,ks)+2.*weight*c(ks)/receptorarea(n)
509    end do
510  end do
511
512  if (mp_measure_time) call mpif_mtime('conccalc',1)
513
514end subroutine conccalc
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG