source: flexpart.git/src/conccalc.f90

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

add SPDX-License-Identifier to all .f90 files

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