source: flexpart.git/src/releaseparticles.f90 @ 75a4ded

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 75a4ded was 75a4ded, checked in by Sabine <sabine.eckhardt@…>, 7 years ago

BUGFIX: ind_r should be treated the same for 1,3 and 4

  • Property mode set to 100644
File size: 17.1 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 releaseparticles(itime)
23  !                              o
24  !*****************************************************************************
25  !                                                                            *
26  !     This subroutine releases particles from the release locations.         *
27  !                                                                            *
28  !     It searches for a "vacant" storage space and assigns all particle      *
29  !     information to that space. A space is vacant either when no particle   *
30  !     is yet assigned to it, or when it's particle is expired and, thus,     *
31  !     the storage space is made available to a new particle.                 *
32  !                                                                            *
33  !     Author: A. Stohl                                                       *
34  !                                                                            *
35  !     29 June 2002                                                           *
36  !                                                                            *
37  !*****************************************************************************
38  !                                                                            *
39  ! Variables:                                                                 *
40  ! itime [s]            current time                                          *
41  ! ireleasestart, ireleaseend          start and end times of all releases    *
42  ! npart(maxpoint)      number of particles to be released in total           *
43  ! numrel               number of particles to be released during this time   *
44  !                      step                                                  *
45  !                                                                            *
46  !*****************************************************************************
47
48  use point_mod
49  use xmass_mod
50  use par_mod
51  use com_mod
52  use random_mod, only: ran1
53
54  implicit none
55
56  !real xaux,yaux,zaux,ran1,rfraction,xmasssave(maxpoint)
57  real :: xaux,yaux,zaux,rfraction
58  real :: topo,rhoaux(2),r,t,rhoout,ddx,ddy,rddx,rddy,p1,p2,p3,p4
59  real :: dz1,dz2,dz,xtn,ytn,xlonav,timecorrect(maxspec),press,pressold
60  real :: presspart,average_timecorrect
61  integer :: itime,numrel,i,j,k,n,ix,jy,ixp,jyp,ipart,minpart,ii
62  integer :: indz,indzp,kz,ngrid
63  integer :: nweeks,ndayofweek,nhour,jjjjmmdd,ihmmss,mm
64  real(kind=dp) :: juldate,julmonday,jul,jullocal,juldiff
65  real,parameter :: eps=nxmax/3.e5,eps2=1.e-6
66
67  integer :: idummy = -7
68  !save idummy,xmasssave
69  !data idummy/-7/,xmasssave/maxpoint*0./
70
71
72
73  ! Determine the actual date and time in Greenwich (i.e., UTC + correction for daylight savings time)
74  !*****************************************************************************
75
76  julmonday=juldate(19000101,0)          ! this is a Monday
77  jul=bdate+real(itime,kind=dp)/86400._dp    ! this is the current day
78  call caldate(jul,jjjjmmdd,ihmmss)
79  mm=(jjjjmmdd-10000*(jjjjmmdd/10000))/100
80  if ((mm.ge.4).and.(mm.le.9)) jul=jul+1._dp/24._dp   ! daylight savings time in summer
81
82
83  ! For every release point, check whether we are in the release time interval
84  !***************************************************************************
85
86  minpart=1
87  do i=1,numpoint
88    if ((itime.ge.ireleasestart(i)).and. &! are we within release interval?
89         (itime.le.ireleaseend(i))) then
90
91  ! Determine the local day and time
92  !*********************************
93
94      xlonav=xlon0+(xpoint2(i)+xpoint1(i))/2.*dx  ! longitude needed to determine local time
95      if (xlonav.lt.-180.) xlonav=xlonav+360.
96      if (xlonav.gt.180.) xlonav=xlonav-360.
97      jullocal=jul+real(xlonav,kind=dp)/360._dp   ! correct approximately for time zone to obtain local time
98
99      juldiff=jullocal-julmonday
100      nweeks=int(juldiff/7._dp)
101      juldiff=juldiff-real(nweeks,kind=dp)*7._dp
102      ndayofweek=int(juldiff)+1              ! this is the current day of week, starting with Monday
103      nhour=nint((juldiff-real(ndayofweek-1,kind=dp))*24._dp)    ! this is the current hour
104      if (nhour.eq.0) then
105        nhour=24
106        ndayofweek=ndayofweek-1
107        if (ndayofweek.eq.0) ndayofweek=7
108      endif
109
110  ! Calculate a species- and time-dependent correction factor, distinguishing between
111  ! area (those with release starting at surface) and point (release starting above surface) sources
112  ! Also, calculate an average time correction factor (species independent)
113  !*****************************************************************************
114      average_timecorrect=0.
115      do k=1,nspec
116        if (zpoint1(i).gt.0.5) then      ! point source
117          timecorrect(k)=point_hour(k,nhour)*point_dow(k,ndayofweek)
118        else                             ! area source
119          timecorrect(k)=area_hour(k,nhour)*area_dow(k,ndayofweek)
120        endif
121        average_timecorrect=average_timecorrect+timecorrect(k)
122      end do
123      average_timecorrect=average_timecorrect/real(nspec)
124
125  ! Determine number of particles to be released this time; at start and at end of release,
126  ! only half the particles are released
127  !*****************************************************************************
128
129      if (ireleasestart(i).ne.ireleaseend(i)) then
130        rfraction=abs(real(npart(i))*real(lsynctime)/ &
131             real(ireleaseend(i)-ireleasestart(i)))
132        if ((itime.eq.ireleasestart(i)).or. &
133             (itime.eq.ireleaseend(i))) rfraction=rfraction/2.
134
135  ! Take the species-average time correction factor in order to scale the
136  ! number of particles released this time
137  !**********************************************************************
138        rfraction=rfraction*average_timecorrect
139
140        rfraction=rfraction+xmasssave(i)  ! number to be released at this time
141        numrel=int(rfraction)
142        xmasssave(i)=rfraction-real(numrel)
143      else
144        numrel=npart(i)
145      endif
146
147      xaux=xpoint2(i)-xpoint1(i)
148      yaux=ypoint2(i)-ypoint1(i)
149      zaux=zpoint2(i)-zpoint1(i)
150      do j=1,numrel                       ! loop over particles to be released this time
151        do ipart=minpart,maxpart          ! search for free storage space
152
153  ! If a free storage space is found, attribute everything to this array element
154  !*****************************************************************************
155
156          if (itra1(ipart).ne.itime) then
157
158  ! Particle coordinates are determined by using a random position within the release volume
159  !*****************************************************************************
160
161  ! Determine horizontal particle position
162  !***************************************
163
164            xtra1(ipart)=xpoint1(i)+ran1(idummy)*xaux
165            if (xglobal) then
166              if (xtra1(ipart).gt.real(nxmin1)) xtra1(ipart)= &
167                   xtra1(ipart)-real(nxmin1)
168              if (xtra1(ipart).lt.0.) xtra1(ipart)= &
169                   xtra1(ipart)+real(nxmin1)
170            endif
171            ytra1(ipart)=ypoint1(i)+ran1(idummy)*yaux
172
173  ! Assign mass to particle: Total mass divided by total number of particles.
174  ! Time variation has partly been taken into account already by a species-average
175  ! correction factor, by which the number of particles released this time has been
176  ! scaled. Adjust the mass per particle by the species-dependent time correction factor
177  ! divided by the species-average one
178  ! for the scavenging calculation the mass needs to be multiplied with rho of the particle layer and
179  ! divided by the sum of rho of all particles.
180  !*****************************************************************************
181            do k=1,nspec
182               xmass1(ipart,k)=xmass(i,k)/real(npart(i)) &
183                    *timecorrect(k)/average_timecorrect
184              if (DRYBKDEP.or.WETBKDEP) then ! if there is no scavenging in wetdepo it will be set to 0
185!              if ( henry(k).gt.0 .or. &
186!                   crain_aero(k).gt.0. .or. csnow_aero(k).gt.0. .or. &
187!                   ccn_aero(k).gt.0. .or. in_aero(k).gt.0. )  then
188                xscav_frac1(ipart,k)=-1.
189               endif
190  ! Assign certain properties to particle
191  !**************************************
192            end do
193            nclass(ipart)=min(int(ran1(idummy)*real(nclassunc))+1, &
194                 nclassunc)
195            numparticlecount=numparticlecount+1
196            if (mquasilag.eq.0) then
197              npoint(ipart)=i
198            else
199              npoint(ipart)=numparticlecount
200            endif
201            idt(ipart)=mintime               ! first time step
202            itra1(ipart)=itime
203            itramem(ipart)=itra1(ipart)
204            itrasplit(ipart)=itra1(ipart)+ldirect*itsplit
205
206
207  ! Determine vertical particle position
208  !*************************************
209
210            ztra1(ipart)=zpoint1(i)+ran1(idummy)*zaux
211  ! Interpolation of topography and density
212  !****************************************
213
214  ! Determine the nest we are in
215  !*****************************
216
217            ngrid=0
218            do k=numbnests,1,-1
219              if ((xtra1(ipart).gt.xln(k)+eps).and. &
220                   (xtra1(ipart).lt.xrn(k)-eps).and. &
221                   (ytra1(ipart).gt.yln(k)+eps).and. &
222                   (ytra1(ipart).lt.yrn(k)-eps)) then
223                ngrid=k
224                goto 43
225              endif
226            end do
22743          continue
228
229  ! Determine (nested) grid coordinates and auxiliary parameters used for interpolation
230  !*****************************************************************************
231
232            if (ngrid.gt.0) then
233              xtn=(xtra1(ipart)-xln(ngrid))*xresoln(ngrid)
234              ytn=(ytra1(ipart)-yln(ngrid))*yresoln(ngrid)
235              ix=int(xtn)
236              jy=int(ytn)
237              ddy=ytn-real(jy)
238              ddx=xtn-real(ix)
239            else
240              ix=int(xtra1(ipart))
241              jy=int(ytra1(ipart))
242              ddy=ytra1(ipart)-real(jy)
243              ddx=xtra1(ipart)-real(ix)
244            endif
245            ixp=ix+1
246            jyp=jy+1
247            rddx=1.-ddx
248            rddy=1.-ddy
249            p1=rddx*rddy
250            p2=ddx*rddy
251            p3=rddx*ddy
252            p4=ddx*ddy
253
254            if (ngrid.gt.0) then
255              topo=p1*oron(ix ,jy ,ngrid) &
256                   + p2*oron(ixp,jy ,ngrid) &
257                   + p3*oron(ix ,jyp,ngrid) &
258                   + p4*oron(ixp,jyp,ngrid)
259            else
260              topo=p1*oro(ix ,jy) &
261                   + p2*oro(ixp,jy) &
262                   + p3*oro(ix ,jyp) &
263                   + p4*oro(ixp,jyp)
264            endif
265
266  ! If starting height is in pressure coordinates, retrieve pressure profile and convert zpart1 to meters
267  !*****************************************************************************
268            if (kindz(i).eq.3) then
269              presspart=ztra1(ipart)
270              do kz=1,nz
271                if (ngrid.gt.0) then
272                  r=p1*rhon(ix ,jy ,kz,2,ngrid) &
273                       +p2*rhon(ixp,jy ,kz,2,ngrid) &
274                       +p3*rhon(ix ,jyp,kz,2,ngrid) &
275                       +p4*rhon(ixp,jyp,kz,2,ngrid)
276                  t=p1*ttn(ix ,jy ,kz,2,ngrid) &
277                       +p2*ttn(ixp,jy ,kz,2,ngrid) &
278                       +p3*ttn(ix ,jyp,kz,2,ngrid) &
279                       +p4*ttn(ixp,jyp,kz,2,ngrid)
280                else
281                  r=p1*rho(ix ,jy ,kz,2) &
282                       +p2*rho(ixp,jy ,kz,2) &
283                       +p3*rho(ix ,jyp,kz,2) &
284                       +p4*rho(ixp,jyp,kz,2)
285                  t=p1*tt(ix ,jy ,kz,2) &
286                       +p2*tt(ixp,jy ,kz,2) &
287                       +p3*tt(ix ,jyp,kz,2) &
288                       +p4*tt(ixp,jyp,kz,2)
289                endif
290                press=r*r_air*t/100.
291                if (kz.eq.1) pressold=press
292
293                if (press.lt.presspart) then
294                  if (kz.eq.1) then
295                    ztra1(ipart)=height(1)/2.
296                  else
297                    dz1=pressold-presspart
298                    dz2=presspart-press
299                    ztra1(ipart)=(height(kz-1)*dz2+height(kz)*dz1) &
300                         /(dz1+dz2)
301                  endif
302                  goto 71
303                endif
304                pressold=press
305              end do
30671            continue
307            endif
308
309  ! If release positions are given in meters above sea level, subtract the
310  ! topography from the starting height
311  !***********************************************************************
312
313            if (kindz(i).eq.2) ztra1(ipart)=ztra1(ipart)-topo
314            if (ztra1(ipart).lt.eps2) ztra1(ipart)=eps2   ! Minimum starting height is eps2
315            if (ztra1(ipart).gt.height(nz)-0.5) ztra1(ipart)= &
316                 height(nz)-0.5 ! Maximum starting height is uppermost level - 0.5 meters
317
318
319
320  ! For special simulations, multiply particle concentration air density;
321  ! Simply take the 2nd field in memory to do this (accurate enough)
322  !***********************************************************************
323  !AF IND_SOURCE switches between different units for concentrations at the source
324  !Af    NOTE that in backward simulations the release of particles takes place at the
325  !Af         receptor and the sampling at the source.
326  !Af          1="mass"
327  !Af          2="mass mixing ratio"
328  !Af IND_RECEPTOR switches between different units for concentrations at the receptor
329  !Af          1="mass"
330  !Af          2="mass mixing ratio"
331
332  !Af switches for the releasefile:
333  !Af IND_REL =  1 : xmass * rho
334  !Af IND_REL =  0 : xmass * 1
335
336  !Af ind_rel is defined in readcommand.f
337
338            if ((ind_rel .eq. 1).or.(ind_rel .eq. 3).or.(ind_rel .eq. 4)) then
339
340  ! Interpolate the air density
341  !****************************
342
343              do ii=2,nz
344                if (height(ii).gt.ztra1(ipart)) then
345                  indz=ii-1
346                  indzp=ii
347                  goto 6
348                endif
349              end do
3506             continue
351
352              dz1=ztra1(ipart)-height(indz)
353              dz2=height(indzp)-ztra1(ipart)
354              dz=1./(dz1+dz2)
355
356              if (ngrid.gt.0) then
357                do n=1,2
358                  rhoaux(n)=p1*rhon(ix ,jy ,indz+n-1,2,ngrid) &
359                       +p2*rhon(ixp,jy ,indz+n-1,2,ngrid) &
360                       +p3*rhon(ix ,jyp,indz+n-1,2,ngrid) &
361                       +p4*rhon(ixp,jyp,indz+n-1,2,ngrid)
362                end do
363              else
364                do n=1,2
365                  rhoaux(n)=p1*rho(ix ,jy ,indz+n-1,2) &
366                       +p2*rho(ixp,jy ,indz+n-1,2) &
367                       +p3*rho(ix ,jyp,indz+n-1,2) &
368                       +p4*rho(ixp,jyp,indz+n-1,2)
369                end do
370              endif
371              rhoout=(dz2*rhoaux(1)+dz1*rhoaux(2))*dz
372              rho_rel(i)=rhoout
373
374
375  ! Multiply "mass" (i.e., mass mixing ratio in forward runs) with density
376  !********************************************************************
377
378              do k=1,nspec
379                xmass1(ipart,k)=xmass1(ipart,k)*rhoout
380              end do
381            endif
382
383            numpart=max(numpart,ipart)
384            goto 34      ! Storage space has been found, stop searching
385          endif
386        end do  ! i=1:numpoint
387        if (ipart.gt.maxpart) goto 996
388
38934      minpart=ipart+1
390      end do ! ipart=minpart,maxpart
391      endif ! j=1,numrel
392  end do
393
394
395  return
396
397996   continue
398  write(*,*) '#####################################################'
399  write(*,*) '#### FLEXPART MODEL SUBROUTINE RELEASEPARTICLES: ####'
400  write(*,*) '####                                             ####'
401  write(*,*) '#### ERROR - TOTAL NUMBER OF PARTICLES REQUIRED  ####'
402  write(*,*) '#### EXCEEDS THE MAXIMUM ALLOWED NUMBER. REDUCE  ####'
403  write(*,*) '#### EITHER NUMBER OF PARTICLES PER RELEASE POINT####'
404  write(*,*) '#### OR REDUCE NUMBER OF RELEASE POINTS.         ####'
405  write(*,*) '#####################################################'
406  stop
407
408end subroutine releaseparticles
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG