source: branches/jerome/src_flexwrf_v3.1/releaseparticles_irreg.f90 @ 16

Last change on this file since 16 was 16, checked in by jebri, 11 years ago

sources for flexwrf v3.1

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