source: flexpart.git/src/boundcond_domainfill.f90

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file was 92fab65, checked in by Ignacio Pisso <ip@…>, 4 years ago

add SPDX-License-Identifier to all .f90 files

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