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