source: flexpart.git/src/boundcond_domainfill.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: 19.3 KB
Line 
1subroutine boundcond_domainfill(itime,loutend)
2  !                                  i      i
3  !*****************************************************************************
4  !                                                                            *
5  ! Particles are created by this subroutine continuously throughout the       *
6  ! simulation at the boundaries of the domain-filling box.                    *
7  ! All particles carry the same amount of mass which alltogether comprises the*
8  ! mass of air within the box, which remains (more or less) constant.         *
9  !                                                                            *
10  !     Author: A. Stohl                                                       *
11  !                                                                            *
12  !     16 October 2002                                                        *
13  !                                                                            *
14  !*****************************************************************************
15  !                                                                            *
16  ! Variables:                                                                 *
17  !                                                                            *
18  ! nx_we(2)       grid indices for western and eastern boundary of domain-    *
19  !                filling trajectory calculations                             *
20  ! ny_sn(2)       grid indices for southern and northern boundary of domain-  *
21  !                filling trajectory calculations                             *
22  !                                                                            *
23  !*****************************************************************************
24
25  use point_mod
26  use par_mod
27  use com_mod
28  use random_mod, only: ran1
29
30  implicit none
31
32  real :: dz,dz1,dz2,dt1,dt2,dtt,ylat,xm,cosfact,accmasst
33  integer :: itime,in,indz,indzp,i,loutend
34  integer :: j,k,ix,jy,m,indzh,indexh,minpart,ipart,mmass
35  integer :: numactiveparticles
36
37  real :: windl(2),rhol(2)
38  real :: windhl(2),rhohl(2)
39  real :: windx,rhox
40  real :: deltaz,boundarea,fluxofmass
41
42  integer :: ixm,ixp,jym,jyp,indzm,mm
43  real :: pvpart,ddx,ddy,rddx,rddy,p1,p2,p3,p4,y1(2),yh1(2)
44
45  integer :: idummy = -11
46
47
48  ! If domain-filling is global, no boundary conditions are needed
49  !***************************************************************
50
51  if (gdomainfill) return
52
53  accmasst=0.
54  numactiveparticles=0
55
56  ! Terminate trajectories that have left the domain, if domain-filling
57  ! trajectory calculation domain is not global
58  !********************************************************************
59
60  do i=1,numpart
61    if (itra1(i).eq.itime) then
62      if ((ytra1(i).gt.real(ny_sn(2))).or. &
63           (ytra1(i).lt.real(ny_sn(1)))) itra1(i)=-999999999
64      if (((.not.xglobal).or.(nx_we(2).ne.(nx-2))).and. &
65           ((xtra1(i).lt.real(nx_we(1))).or. &
66           (xtra1(i).gt.real(nx_we(2))))) itra1(i)=-999999999
67    endif
68    if (itra1(i).ne.-999999999) numactiveparticles= &
69         numactiveparticles+1
70  end do
71
72
73  ! Determine auxiliary variables for time interpolation
74  !*****************************************************
75
76  dt1=real(itime-memtime(1))
77  dt2=real(memtime(2)-itime)
78  dtt=1./(dt1+dt2)
79
80  ! Initialize auxiliary variable used to search for vacant storage space
81  !**********************************************************************
82
83  minpart=1
84
85  !***************************************
86  ! Western and eastern boundary condition
87  !***************************************
88
89  ! Loop from south to north
90  !*************************
91
92  do jy=ny_sn(1),ny_sn(2)
93
94  ! Loop over western (index 1) and eastern (index 2) boundary
95  !***********************************************************
96
97    do k=1,2
98
99  ! Loop over all release locations in a column
100  !********************************************
101
102      do j=1,numcolumn_we(k,jy)
103
104  ! Determine, for each release location, the area of the corresponding boundary
105  !*****************************************************************************
106
107        if (j.eq.1) then
108          deltaz=(zcolumn_we(k,jy,2)+zcolumn_we(k,jy,1))/2.
109        else if (j.eq.numcolumn_we(k,jy)) then
110  !        deltaz=height(nz)-(zcolumn_we(k,jy,j-1)+
111  !    +        zcolumn_we(k,jy,j))/2.
112  ! In order to avoid taking a very high column for very many particles,
113  ! use the deltaz from one particle below instead
114          deltaz=(zcolumn_we(k,jy,j)-zcolumn_we(k,jy,j-2))/2.
115        else
116          deltaz=(zcolumn_we(k,jy,j+1)-zcolumn_we(k,jy,j-1))/2.
117        endif
118        if ((jy.eq.ny_sn(1)).or.(jy.eq.ny_sn(2))) then
119          boundarea=deltaz*111198.5/2.*dy
120        else
121          boundarea=deltaz*111198.5*dy
122        endif
123
124
125  ! Interpolate the wind velocity and density to the release location
126  !******************************************************************
127
128  ! Determine the model level below the release position
129  !*****************************************************
130
131        do i=2,nz
132          if (height(i).gt.zcolumn_we(k,jy,j)) then
133            indz=i-1
134            indzp=i
135            goto 6
136          endif
137        end do
1386       continue
139
140  ! Vertical distance to the level below and above current position
141  !****************************************************************
142
143        dz1=zcolumn_we(k,jy,j)-height(indz)
144        dz2=height(indzp)-zcolumn_we(k,jy,j)
145        dz=1./(dz1+dz2)
146
147  ! Vertical and temporal interpolation
148  !************************************
149
150        do m=1,2
151          indexh=memind(m)
152          do in=1,2
153            indzh=indz+in-1
154            windl(in)=uu(nx_we(k),jy,indzh,indexh)
155            rhol(in)=rho(nx_we(k),jy,indzh,indexh)
156          end do
157
158          windhl(m)=(dz2*windl(1)+dz1*windl(2))*dz
159          rhohl(m)=(dz2*rhol(1)+dz1*rhol(2))*dz
160        end do
161
162        windx=(windhl(1)*dt2+windhl(2)*dt1)*dtt
163        rhox=(rhohl(1)*dt2+rhohl(2)*dt1)*dtt
164
165  ! Calculate mass flux
166  !********************
167
168        fluxofmass=windx*rhox*boundarea*real(lsynctime)
169
170
171  ! If the mass flux is directed into the domain, add it to previous mass fluxes;
172  ! if it is out of the domain, set accumulated mass flux to zero
173  !******************************************************************************
174
175        if (k.eq.1) then
176          if (fluxofmass.ge.0.) then
177            acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)+fluxofmass
178          else
179            acc_mass_we(k,jy,j)=0.
180          endif
181        else
182          if (fluxofmass.le.0.) then
183            acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)+abs(fluxofmass)
184          else
185            acc_mass_we(k,jy,j)=0.
186          endif
187        endif
188        accmasst=accmasst+acc_mass_we(k,jy,j)
189
190  ! If the accumulated mass exceeds half the mass that each particle shall carry,
191  ! one (or more) particle(s) is (are) released and the accumulated mass is
192  ! reduced by the mass of this (these) particle(s)
193  !******************************************************************************
194
195        if (acc_mass_we(k,jy,j).ge.xmassperparticle/2.) then
196          mmass=int((acc_mass_we(k,jy,j)+xmassperparticle/2.)/ &
197               xmassperparticle)
198          acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)- &
199               real(mmass)*xmassperparticle
200        else
201          mmass=0
202        endif
203
204        do m=1,mmass
205          do ipart=minpart,maxpart
206
207  ! If a vacant storage space is found, attribute everything to this array element
208  !*****************************************************************************
209
210            if (itra1(ipart).ne.itime) then
211
212  ! Assign particle positions
213  !**************************
214
215              xtra1(ipart)=real(nx_we(k))
216              if (jy.eq.ny_sn(1)) then
217                ytra1(ipart)=real(jy)+0.5*ran1(idummy)
218              else if (jy.eq.ny_sn(2)) then
219                ytra1(ipart)=real(jy)-0.5*ran1(idummy)
220              else
221                ytra1(ipart)=real(jy)+(ran1(idummy)-.5)
222              endif
223              if (j.eq.1) then
224                ztra1(ipart)=zcolumn_we(k,jy,1)+(zcolumn_we(k,jy,2)- &
225                     zcolumn_we(k,jy,1))/4.
226              else if (j.eq.numcolumn_we(k,jy)) then
227                ztra1(ipart)=(2.*zcolumn_we(k,jy,j)+ &
228                     zcolumn_we(k,jy,j-1)+height(nz))/4.
229              else
230                ztra1(ipart)=zcolumn_we(k,jy,j-1)+ran1(idummy)* &
231                     (zcolumn_we(k,jy,j+1)-zcolumn_we(k,jy,j-1))
232              endif
233
234  ! Interpolate PV to the particle position
235  !****************************************
236              ixm=int(xtra1(ipart))
237              jym=int(ytra1(ipart))
238              ixp=ixm+1
239              jyp=jym+1
240              ddx=xtra1(ipart)-real(ixm)
241              ddy=ytra1(ipart)-real(jym)
242              rddx=1.-ddx
243              rddy=1.-ddy
244              p1=rddx*rddy
245              p2=ddx*rddy
246              p3=rddx*ddy
247              p4=ddx*ddy
248              do i=2,nz
249                if (height(i).gt.ztra1(ipart)) then
250                  indzm=i-1
251                  indzp=i
252                  goto 26
253                endif
254              end do
25526            continue
256              dz1=ztra1(ipart)-height(indzm)
257              dz2=height(indzp)-ztra1(ipart)
258              dz=1./(dz1+dz2)
259              do mm=1,2
260                indexh=memind(mm)
261                do in=1,2
262                  indzh=indzm+in-1
263                  y1(in)=p1*pv(ixm,jym,indzh,indexh) &
264                       +p2*pv(ixp,jym,indzh,indexh) &
265                       +p3*pv(ixm,jyp,indzh,indexh) &
266                       +p4*pv(ixp,jyp,indzh,indexh)
267                end do
268                yh1(mm)=(dz2*y1(1)+dz1*y1(2))*dz
269              end do
270              pvpart=(yh1(1)*dt2+yh1(2)*dt1)*dtt
271              ylat=ylat0+ytra1(ipart)*dy
272              if (ylat.lt.0.) pvpart=-1.*pvpart
273
274
275  ! For domain-filling option 2 (stratospheric O3), do the rest only in the stratosphere
276  !*****************************************************************************
277
278              if (((ztra1(ipart).gt.3000.).and. &
279                   (pvpart.gt.pvcrit)).or.(mdomainfill.eq.1)) then
280                nclass(ipart)=min(int(ran1(idummy)* &
281                     real(nclassunc))+1,nclassunc)
282                numactiveparticles=numactiveparticles+1
283                numparticlecount=numparticlecount+1
284                npoint(ipart)=numparticlecount
285                idt(ipart)=mintime
286                itra1(ipart)=itime
287                itramem(ipart)=itra1(ipart)
288                itrasplit(ipart)=itra1(ipart)+ldirect*itsplit
289                xmass1(ipart,1)=xmassperparticle
290                if (mdomainfill.eq.2) xmass1(ipart,1)= &
291                     xmass1(ipart,1)*pvpart*48./29.*ozonescale/10.**9
292              else
293                goto 71
294              endif
295
296
297  ! Increase numpart, if necessary
298  !*******************************
299
300              numpart=max(numpart,ipart)
301              goto 73      ! Storage space has been found, stop searching
302            endif
303          end do
304          if (ipart.gt.maxpart) &
305               stop 'boundcond_domainfill.f: too many particles required'
30673        minpart=ipart+1
30771        continue
308        end do
309
310
311      end do
312    end do
313  end do
314
315
316  !*****************************************
317  ! Southern and northern boundary condition
318  !*****************************************
319
320  ! Loop from west to east
321  !***********************
322
323  do ix=nx_we(1),nx_we(2)
324
325  ! Loop over southern (index 1) and northern (index 2) boundary
326  !*************************************************************
327
328    do k=1,2
329      ylat=ylat0+real(ny_sn(k))*dy
330      cosfact=cos(ylat*pi180)
331
332  ! Loop over all release locations in a column
333  !********************************************
334
335      do j=1,numcolumn_sn(k,ix)
336
337  ! Determine, for each release location, the area of the corresponding boundary
338  !*****************************************************************************
339
340        if (j.eq.1) then
341          deltaz=(zcolumn_sn(k,ix,2)+zcolumn_sn(k,ix,1))/2.
342        else if (j.eq.numcolumn_sn(k,ix)) then
343  !        deltaz=height(nz)-(zcolumn_sn(k,ix,j-1)+
344  !    +        zcolumn_sn(k,ix,j))/2.
345  ! In order to avoid taking a very high column for very many particles,
346  ! use the deltaz from one particle below instead
347          deltaz=(zcolumn_sn(k,ix,j)-zcolumn_sn(k,ix,j-2))/2.
348        else
349          deltaz=(zcolumn_sn(k,ix,j+1)-zcolumn_sn(k,ix,j-1))/2.
350        endif
351        if ((ix.eq.nx_we(1)).or.(ix.eq.nx_we(2))) then
352          boundarea=deltaz*111198.5/2.*cosfact*dx
353        else
354          boundarea=deltaz*111198.5*cosfact*dx
355        endif
356
357
358  ! Interpolate the wind velocity and density to the release location
359  !******************************************************************
360
361  ! Determine the model level below the release position
362  !*****************************************************
363
364        do i=2,nz
365          if (height(i).gt.zcolumn_sn(k,ix,j)) then
366            indz=i-1
367            indzp=i
368            goto 16
369          endif
370        end do
37116      continue
372
373  ! Vertical distance to the level below and above current position
374  !****************************************************************
375
376        dz1=zcolumn_sn(k,ix,j)-height(indz)
377        dz2=height(indzp)-zcolumn_sn(k,ix,j)
378        dz=1./(dz1+dz2)
379
380  ! Vertical and temporal interpolation
381  !************************************
382
383        do m=1,2
384          indexh=memind(m)
385          do in=1,2
386            indzh=indz+in-1
387            windl(in)=vv(ix,ny_sn(k),indzh,indexh)
388            rhol(in)=rho(ix,ny_sn(k),indzh,indexh)
389          end do
390
391          windhl(m)=(dz2*windl(1)+dz1*windl(2))*dz
392          rhohl(m)=(dz2*rhol(1)+dz1*rhol(2))*dz
393        end do
394
395        windx=(windhl(1)*dt2+windhl(2)*dt1)*dtt
396        rhox=(rhohl(1)*dt2+rhohl(2)*dt1)*dtt
397
398  ! Calculate mass flux
399  !********************
400
401        fluxofmass=windx*rhox*boundarea*real(lsynctime)
402
403  ! If the mass flux is directed into the domain, add it to previous mass fluxes;
404  ! if it is out of the domain, set accumulated mass flux to zero
405  !******************************************************************************
406
407        if (k.eq.1) then
408          if (fluxofmass.ge.0.) then
409            acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)+fluxofmass
410          else
411            acc_mass_sn(k,ix,j)=0.
412          endif
413        else
414          if (fluxofmass.le.0.) then
415            acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)+abs(fluxofmass)
416          else
417            acc_mass_sn(k,ix,j)=0.
418          endif
419        endif
420        accmasst=accmasst+acc_mass_sn(k,ix,j)
421
422  ! If the accumulated mass exceeds half the mass that each particle shall carry,
423  ! one (or more) particle(s) is (are) released and the accumulated mass is
424  ! reduced by the mass of this (these) particle(s)
425  !******************************************************************************
426
427        if (acc_mass_sn(k,ix,j).ge.xmassperparticle/2.) then
428          mmass=int((acc_mass_sn(k,ix,j)+xmassperparticle/2.)/ &
429               xmassperparticle)
430          acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)- &
431               real(mmass)*xmassperparticle
432        else
433          mmass=0
434        endif
435
436        do m=1,mmass
437          do ipart=minpart,maxpart
438
439  ! If a vacant storage space is found, attribute everything to this array element
440  !*****************************************************************************
441
442            if (itra1(ipart).ne.itime) then
443
444  ! Assign particle positions
445  !**************************
446
447              ytra1(ipart)=real(ny_sn(k))
448              if (ix.eq.nx_we(1)) then
449                xtra1(ipart)=real(ix)+0.5*ran1(idummy)
450              else if (ix.eq.nx_we(2)) then
451                xtra1(ipart)=real(ix)-0.5*ran1(idummy)
452              else
453                xtra1(ipart)=real(ix)+(ran1(idummy)-.5)
454              endif
455              if (j.eq.1) then
456                ztra1(ipart)=zcolumn_sn(k,ix,1)+(zcolumn_sn(k,ix,2)- &
457                     zcolumn_sn(k,ix,1))/4.
458              else if (j.eq.numcolumn_sn(k,ix)) then
459                ztra1(ipart)=(2.*zcolumn_sn(k,ix,j)+ &
460                     zcolumn_sn(k,ix,j-1)+height(nz))/4.
461              else
462                ztra1(ipart)=zcolumn_sn(k,ix,j-1)+ran1(idummy)* &
463                     (zcolumn_sn(k,ix,j+1)-zcolumn_sn(k,ix,j-1))
464              endif
465
466
467  ! Interpolate PV to the particle position
468  !****************************************
469              ixm=int(xtra1(ipart))
470              jym=int(ytra1(ipart))
471              ixp=ixm+1
472              jyp=jym+1
473              ddx=xtra1(ipart)-real(ixm)
474              ddy=ytra1(ipart)-real(jym)
475              rddx=1.-ddx
476              rddy=1.-ddy
477              p1=rddx*rddy
478              p2=ddx*rddy
479              p3=rddx*ddy
480              p4=ddx*ddy
481              do i=2,nz
482                if (height(i).gt.ztra1(ipart)) then
483                  indzm=i-1
484                  indzp=i
485                  goto 126
486                endif
487              end do
488126           continue
489              dz1=ztra1(ipart)-height(indzm)
490              dz2=height(indzp)-ztra1(ipart)
491              dz=1./(dz1+dz2)
492              do mm=1,2
493                indexh=memind(mm)
494                do in=1,2
495                  indzh=indzm+in-1
496                  y1(in)=p1*pv(ixm,jym,indzh,indexh) &
497                       +p2*pv(ixp,jym,indzh,indexh) &
498                       +p3*pv(ixm,jyp,indzh,indexh) &
499                       +p4*pv(ixp,jyp,indzh,indexh)
500                end do
501                yh1(mm)=(dz2*y1(1)+dz1*y1(2))*dz
502              end do
503              pvpart=(yh1(1)*dt2+yh1(2)*dt1)*dtt
504              if (ylat.lt.0.) pvpart=-1.*pvpart
505
506
507  ! For domain-filling option 2 (stratospheric O3), do the rest only in the stratosphere
508  !*****************************************************************************
509
510              if (((ztra1(ipart).gt.3000.).and. &
511                   (pvpart.gt.pvcrit)).or.(mdomainfill.eq.1)) then
512                nclass(ipart)=min(int(ran1(idummy)* &
513                     real(nclassunc))+1,nclassunc)
514                numactiveparticles=numactiveparticles+1
515                numparticlecount=numparticlecount+1
516                npoint(ipart)=numparticlecount
517                idt(ipart)=mintime
518                itra1(ipart)=itime
519                itramem(ipart)=itra1(ipart)
520                itrasplit(ipart)=itra1(ipart)+ldirect*itsplit
521                xmass1(ipart,1)=xmassperparticle
522                if (mdomainfill.eq.2) xmass1(ipart,1)= &
523                     xmass1(ipart,1)*pvpart*48./29.*ozonescale/10.**9
524              else
525                goto 171
526              endif
527
528
529  ! Increase numpart, if necessary
530  !*******************************
531              numpart=max(numpart,ipart)
532              goto 173      ! Storage space has been found, stop searching
533            endif
534          end do
535          if (ipart.gt.maxpart) &
536               stop 'boundcond_domainfill.f: too many particles required'
537173       minpart=ipart+1
538171       continue
539        end do
540
541
542      end do
543    end do
544  end do
545
546
547  xm=0.
548  do i=1,numpart
549    if (itra1(i).eq.itime) xm=xm+xmass1(i,1)
550  end do
551
552  !write(*,*) itime,numactiveparticles,numparticlecount,numpart,
553  !    +xm,accmasst,xm+accmasst
554
555
556  ! If particles shall be dumped, then accumulated masses at the domain boundaries
557  ! must be dumped, too, to be used for later runs
558  !*****************************************************************************
559
560  if ((ipout.gt.0).and.(itime.eq.loutend)) then
561    open(unitboundcond,file=path(2)(1:length(2))//'boundcond.bin', &
562         form='unformatted')
563    write(unitboundcond) numcolumn_we,numcolumn_sn, &
564         zcolumn_we,zcolumn_sn,acc_mass_we,acc_mass_sn
565    close(unitboundcond)
566  endif
567
568end subroutine boundcond_domainfill
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG