releaseparticles_mpi.f90 17.8 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
!**********************************************************************
! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010         *
! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa,             *
! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann   *
!                                                                     *
! This file is part of FLEXPART.                                      *
!                                                                     *
! FLEXPART is free software: you can redistribute it and/or modify    *
! it under the terms of the GNU General Public License as published by*
! the Free Software Foundation, either version 3 of the License, or   *
! (at your option) any later version.                                 *
!                                                                     *
! FLEXPART is distributed in the hope that it will be useful,         *
! but WITHOUT ANY WARRANTY; without even the implied warranty of      *
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
! GNU General Public License for more details.                        *
!                                                                     *
! You should have received a copy of the GNU General Public License   *
! along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
!**********************************************************************

subroutine releaseparticles(itime)
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
!                              o
!*****************************************************************************
!                                                                            *
!     This subroutine releases particles from the release locations.         *
!                                                                            *
!     It searches for a "vacant" storage space and assigns all particle      *
!     information to that space. A space is vacant either when no particle   *
!     is yet assigned to it, or when it's particle is expired and, thus,     *
!     the storage space is made available to a new particle.                 *
!                                                                            *
!     Author: A. Stohl                                                       *
!                                                                            *
!     29 June 2002                                                           *
!                                                                            *
!   CHANGES                                                                  *
!     12/2014 eso: MPI version                                               *
!                                                                            *
!*****************************************************************************
!                                                                            *
! Variables:                                                                 *
! itime [s]            current time                                          *
! ireleasestart, ireleaseend          start and end times of all releases    *
! npart(maxpoint)      number of particles to be released in total           *
! numrel               number of particles to be released during this time   *
!                      step                                                  *
! addone               extra particle assigned to MPI process if             *
!                      mod(npart(i),mp_partgroup_np) .ne. 0)                 *
!*****************************************************************************
51
52
53
54
55
56

  use point_mod
  use xmass_mod
  use par_mod
  use com_mod
  use random_mod, only: ran1
57
  use mpi_mod, only: mp_partid, maxpart_mpi, mp_partgroup_np, mp_seed, mpif_mpi_barrier
58
59
60

  implicit none

61
!real xaux,yaux,zaux,ran1,rfraction,xmasssave(maxpoint)
62
63
64
65
66
67
68
69
70
  real :: xaux,yaux,zaux,rfraction
  real :: topo,rhoaux(2),r,t,rhoout,ddx,ddy,rddx,rddy,p1,p2,p3,p4
  real :: dz1,dz2,dz,xtn,ytn,xlonav,timecorrect(maxspec),press,pressold
  real :: presspart,average_timecorrect
  integer :: itime,numrel,i,j,k,n,ix,jy,ixp,jyp,ipart,minpart,ii,addone
  integer :: indz,indzp,kz,ngrid
  integer :: nweeks,ndayofweek,nhour,jjjjmmdd,ihmmss,mm
  real(kind=dp) :: juldate,julmonday,jul,jullocal,juldiff
  real,parameter :: eps=nxmax/3.e5,eps2=1.e-6
Espen Sollum's avatar
Espen Sollum committed
71
72
  integer :: mind2
! mind2        eso: pointer to 2nd windfield in memory
73
74

  integer :: idummy = -7
75
76
!save idummy,xmasssave
!data idummy/-7/,xmasssave/maxpoint*0./
77
78
79

  logical :: first_call=.true.

80
81
! Use different seed for each process.
!****************************************************************************
82
83
84
85
86
  if (first_call) then
    idummy=idummy+mp_seed
    first_call=.false.
  end if

87
  mind2=memind(2)
88

89
90
! Determine the actual date and time in Greenwich (i.e., UTC + correction for daylight savings time)
!*****************************************************************************
91
92
93
94
95
96
97
98

  julmonday=juldate(19000101,0)          ! this is a Monday
  jul=bdate+real(itime,kind=dp)/86400._dp    ! this is the current day
  call caldate(jul,jjjjmmdd,ihmmss)
  mm=(jjjjmmdd-10000*(jjjjmmdd/10000))/100
  if ((mm.ge.4).and.(mm.le.9)) jul=jul+1._dp/24._dp   ! daylight savings time in summer


99
100
! For every release point, check whether we are in the release time interval
!***************************************************************************
101
102
103
104
105
106

  minpart=1
  do i=1,numpoint
    if ((itime.ge.ireleasestart(i)).and. &! are we within release interval?
         (itime.le.ireleaseend(i))) then

107
108
! Determine the local day and time
!*********************************
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125

      xlonav=xlon0+(xpoint2(i)+xpoint1(i))/2.*dx  ! longitude needed to determine local time
      if (xlonav.lt.-180.) xlonav=xlonav+360.
      if (xlonav.gt.180.) xlonav=xlonav-360.
      jullocal=jul+real(xlonav,kind=dp)/360._dp   ! correct approximately for time zone to obtain local time

      juldiff=jullocal-julmonday
      nweeks=int(juldiff/7._dp)
      juldiff=juldiff-real(nweeks,kind=dp)*7._dp
      ndayofweek=int(juldiff)+1              ! this is the current day of week, starting with Monday
      nhour=nint((juldiff-real(ndayofweek-1,kind=dp))*24._dp)    ! this is the current hour
      if (nhour.eq.0) then
        nhour=24
        ndayofweek=ndayofweek-1
        if (ndayofweek.eq.0) ndayofweek=7
      endif

126
127
128
129
! Calculate a species- and time-dependent correction factor, distinguishing between
! area (those with release starting at surface) and point (release starting above surface) sources
! Also, calculate an average time correction factor (species independent)
!*****************************************************************************
130
131
132
133
134
135
136
137
138
139
140
      average_timecorrect=0.
      do k=1,nspec
        if (zpoint1(i).gt.0.5) then      ! point source
          timecorrect(k)=point_hour(k,nhour)*point_dow(k,ndayofweek)
        else                             ! area source
          timecorrect(k)=area_hour(k,nhour)*area_dow(k,ndayofweek)
        endif
        average_timecorrect=average_timecorrect+timecorrect(k)
      end do
      average_timecorrect=average_timecorrect/real(nspec)

141
142
143
! Determine number of particles to be released this time; at start and at end of release,
! only half the particles are released
!*****************************************************************************
144
145
146
147
148
149
150

      if (ireleasestart(i).ne.ireleaseend(i)) then
        rfraction=abs(real(npart(i))*real(lsynctime)/ &
             real(ireleaseend(i)-ireleasestart(i)))
        if ((itime.eq.ireleasestart(i)).or. &
             (itime.eq.ireleaseend(i))) rfraction=rfraction/2.

151
152
153
154
! Take the species-average time correction factor in order to scale the
! number of particles released this time
! Also scale by number of MPI processes
!**********************************************************************
155

156
157
158
159
        rfraction=rfraction*average_timecorrect

        rfraction=rfraction+xmasssave(i)  ! number to be released at this time

160
! number to be released for one process
161
        if (mp_partid.lt.mod(int(rfraction),mp_partgroup_np)) then
162
163
164
165
166
          addone=1
        else
          addone=0
        end if

167
168
169
        numrel=int(rfraction/mp_partgroup_np) + addone

        xmasssave(i)=rfraction-int(rfraction)
170
171

      else
172
173
! All particles are released in this time interval
! ************************************************
174
175
176
177
178
179
180
181
        if (mp_partid.lt.mod(npart(i),mp_partgroup_np)) then
          addone=1
        else
          addone=0
        end if

        numrel=npart(i)/mp_partgroup_np+addone
      endif
182

183
184
185
186
187
188
      xaux=xpoint2(i)-xpoint1(i)
      yaux=ypoint2(i)-ypoint1(i)
      zaux=zpoint2(i)-zpoint1(i)
      do j=1,numrel                       ! loop over particles to be released this time
        do ipart=minpart,maxpart_mpi     ! search for free storage space

189
190
! If a free storage space is found, attribute everything to this array element
!*****************************************************************************
191
192
193

          if (itra1(ipart).ne.itime) then

194
195
! Particle coordinates are determined by using a random position within the release volume
!*****************************************************************************
196

197
198
! Determine horizontal particle position
!***************************************
199
200
201
202
203
204
205
206
207
208

            xtra1(ipart)=xpoint1(i)+ran1(idummy)*xaux
            if (xglobal) then
              if (xtra1(ipart).gt.real(nxmin1)) xtra1(ipart)= &
                   xtra1(ipart)-real(nxmin1)
              if (xtra1(ipart).lt.0.) xtra1(ipart)= &
                   xtra1(ipart)+real(nxmin1)
            endif
            ytra1(ipart)=ypoint1(i)+ran1(idummy)*yaux

209
210
211
212
213
214
! Assign mass to particle: Total mass divided by total number of particles.
! Time variation has partly been taken into account already by a species-average
! correction factor, by which the number of particles released this time has been
! scaled. Adjust the mass per particle by the species-dependent time correction factor
! divided by the species-average one
!*****************************************************************************
215
            do k=1,nspec
216
217
218
219
220
              xmass1(ipart,k)=xmass(i,k)/real(npart(i)) &
                   *timecorrect(k)/average_timecorrect
!             write (*,*) 'xmass1: ',xmass1(ipart,k),ipart,k
! Assign certain properties to particle
!**************************************
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
            end do
            nclass(ipart)=min(int(ran1(idummy)*real(nclassunc))+1, &
                 nclassunc)
            numparticlecount=numparticlecount+1
            if (mquasilag.eq.0) then
              npoint(ipart)=i
            else
              npoint(ipart)=numparticlecount
            endif
            idt(ipart)=mintime               ! first time step
            itra1(ipart)=itime
            itramem(ipart)=itra1(ipart)
            itrasplit(ipart)=itra1(ipart)+ldirect*itsplit


236
237
! Determine vertical particle position
!*************************************
238
239
240

            ztra1(ipart)=zpoint1(i)+ran1(idummy)*zaux

241
242
! Interpolation of topography and density
!****************************************
243

244
245
! Determine the nest we are in
!*****************************
246
247
248
249
250
251
252
253
254
255
256
257
258

            ngrid=0
            do k=numbnests,1,-1
              if ((xtra1(ipart).gt.xln(k)+eps).and. &
                   (xtra1(ipart).lt.xrn(k)-eps).and. &
                   (ytra1(ipart).gt.yln(k)+eps).and. &
                   (ytra1(ipart).lt.yrn(k)-eps)) then
                ngrid=k
                goto 43
              endif
            end do
43          continue

259
260
! Determine (nested) grid coordinates and auxiliary parameters used for interpolation
!*****************************************************************************
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295

            if (ngrid.gt.0) then
              xtn=(xtra1(ipart)-xln(ngrid))*xresoln(ngrid)
              ytn=(ytra1(ipart)-yln(ngrid))*yresoln(ngrid)
              ix=int(xtn)
              jy=int(ytn)
              ddy=ytn-real(jy)
              ddx=xtn-real(ix)
            else
              ix=int(xtra1(ipart))
              jy=int(ytra1(ipart))
              ddy=ytra1(ipart)-real(jy)
              ddx=xtra1(ipart)-real(ix)
            endif
            ixp=ix+1
            jyp=jy+1
            rddx=1.-ddx
            rddy=1.-ddy
            p1=rddx*rddy
            p2=ddx*rddy
            p3=rddx*ddy
            p4=ddx*ddy

            if (ngrid.gt.0) then
              topo=p1*oron(ix ,jy ,ngrid) &
                   + p2*oron(ixp,jy ,ngrid) &
                   + p3*oron(ix ,jyp,ngrid) &
                   + p4*oron(ixp,jyp,ngrid)
            else
              topo=p1*oro(ix ,jy) &
                   + p2*oro(ixp,jy) &
                   + p3*oro(ix ,jyp) &
                   + p4*oro(ixp,jyp)
            endif

296
297
! If starting height is in pressure coordinates, retrieve pressure profile and convert zpart1 to meters
!*****************************************************************************
298
299
300
301
            if (kindz(i).eq.3) then
              presspart=ztra1(ipart)
              do kz=1,nz
                if (ngrid.gt.0) then
Espen Sollum's avatar
Espen Sollum committed
302
303
304
305
306
307
308
309
                  r=p1*rhon(ix ,jy ,kz,mind2,ngrid) &
                       +p2*rhon(ixp,jy ,kz,mind2,ngrid) &
                       +p3*rhon(ix ,jyp,kz,mind2,ngrid) &
                       +p4*rhon(ixp,jyp,kz,mind2,ngrid)
                  t=p1*ttn(ix ,jy ,kz,mind2,ngrid) &
                       +p2*ttn(ixp,jy ,kz,mind2,ngrid) &
                       +p3*ttn(ix ,jyp,kz,mind2,ngrid) &
                       +p4*ttn(ixp,jyp,kz,mind2,ngrid)
310
                else
Espen Sollum's avatar
Espen Sollum committed
311
312
313
314
315
316
317
318
                  r=p1*rho(ix ,jy ,kz,mind2) &
                       +p2*rho(ixp,jy ,kz,mind2) &
                       +p3*rho(ix ,jyp,kz,mind2) &
                       +p4*rho(ixp,jyp,kz,mind2)
                  t=p1*tt(ix ,jy ,kz,mind2) &
                       +p2*tt(ixp,jy ,kz,mind2) &
                       +p3*tt(ix ,jyp,kz,mind2) &
                       +p4*tt(ixp,jyp,kz,mind2)
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
                endif
                press=r*r_air*t/100.
                if (kz.eq.1) pressold=press

                if (press.lt.presspart) then
                  if (kz.eq.1) then
                    ztra1(ipart)=height(1)/2.
                  else
                    dz1=pressold-presspart
                    dz2=presspart-press
                    ztra1(ipart)=(height(kz-1)*dz2+height(kz)*dz1) &
                         /(dz1+dz2)
                  endif
                  goto 71
                endif
                pressold=press
              end do
71            continue
            endif

339
340
341
! If release positions are given in meters above sea level, subtract the
! topography from the starting height
!***********************************************************************
342
343
344
345
346
347
348
349

            if (kindz(i).eq.2) ztra1(ipart)=ztra1(ipart)-topo
            if (ztra1(ipart).lt.eps2) ztra1(ipart)=eps2   ! Minimum starting height is eps2
            if (ztra1(ipart).gt.height(nz)-0.5) ztra1(ipart)= &
                 height(nz)-0.5 ! Maximum starting height is uppermost level - 0.5 meters



350
351
352
353
354
355
356
357
358
359
360
! For special simulations, multiply particle concentration air density;
! Simply take the 2nd field in memory to do this (accurate enough)
!***********************************************************************
!AF IND_SOURCE switches between different units for concentrations at the source
!Af    NOTE that in backward simulations the release of particles takes place at the
!Af         receptor and the sampling at the source.
!Af          1="mass"
!Af          2="mass mixing ratio"
!Af IND_RECEPTOR switches between different units for concentrations at the receptor
!Af          1="mass"
!Af          2="mass mixing ratio"
361

362
363
364
!Af switches for the releasefile:
!Af IND_REL =  1 : xmass * rho
!Af IND_REL =  0 : xmass * 1
365

366
!Af ind_rel is defined in readcommand.f
367
368
369

            if (ind_rel .eq. 1) then

370
371
! Interpolate the air density
!****************************
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387

              do ii=2,nz
                if (height(ii).gt.ztra1(ipart)) then
                  indz=ii-1
                  indzp=ii
                  goto 6
                endif
              end do
6             continue

              dz1=ztra1(ipart)-height(indz)
              dz2=height(indzp)-ztra1(ipart)
              dz=1./(dz1+dz2)

              if (ngrid.gt.0) then
                do n=1,2
Espen Sollum's avatar
Espen Sollum committed
388
389
390
391
                  rhoaux(n)=p1*rhon(ix ,jy ,indz+n-1,mind2,ngrid) &
                       +p2*rhon(ixp,jy ,indz+n-1,mind2,ngrid) &
                       +p3*rhon(ix ,jyp,indz+n-1,mind2,ngrid) &
                       +p4*rhon(ixp,jyp,indz+n-1,mind2,ngrid)
392
393
394
                end do
              else
                do n=1,2
Espen Sollum's avatar
Espen Sollum committed
395
396
397
398
                  rhoaux(n)=p1*rho(ix ,jy ,indz+n-1,mind2) &
                       +p2*rho(ixp,jy ,indz+n-1,mind2) &
                       +p3*rho(ix ,jyp,indz+n-1,mind2) &
                       +p4*rho(ixp,jyp,indz+n-1,mind2)
399
400
401
402
403
404
                end do
              endif
              rhoout=(dz2*rhoaux(1)+dz1*rhoaux(2))*dz
              rho_rel(i)=rhoout


405
406
! Multiply "mass" (i.e., mass mixing ratio in forward runs) with density
!********************************************************************
407
408
409
410
411
412
413
414
415
416
417

              do k=1,nspec
                xmass1(ipart,k)=xmass1(ipart,k)*rhoout
              end do
            endif


            numpart=max(numpart,ipart)
            goto 34      ! Storage space has been found, stop searching
          endif
        end do
418
419
! ESO TODO: Here we could use dynamic allocation and increase array sizes
        if (ipart.gt.maxpart_mpi) goto 996
420
421
422

34      minpart=ipart+1
      end do
423
    endif
424
425
426
427
428
  end do


  return

429
996 continue
430
431
432
433
434
435
436
437
438
439
440
  write(*,*) '#####################################################'
  write(*,*) '#### FLEXPART MODEL SUBROUTINE RELEASEPARTICLES: ####'
  write(*,*) '####                                             ####'
  write(*,*) '#### ERROR - TOTAL NUMBER OF PARTICLES REQUIRED  ####'
  write(*,*) '#### EXCEEDS THE MAXIMUM ALLOWED NUMBER. REDUCE  ####'
  write(*,*) '#### EITHER NUMBER OF PARTICLES PER RELEASE POINT####'
  write(*,*) '#### OR REDUCE NUMBER OF RELEASE POINTS.         ####'
  write(*,*) '#####################################################'
  stop

end subroutine releaseparticles