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