source: flexpart.git/src/plumetraj.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: 7.5 KB
Line 
1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
3
4subroutine plumetraj(itime)
5  !                       i
6  !*****************************************************************************
7  !                                                                            *
8  ! Determines a plume centroid trajectory for each release site, and manages  *
9  ! clustering of particle locations. Certain parameters (average PV,          *
10  ! tropopause height, etc., are provided along the plume trajectories.        *
11  ! At the end, output is written to file 'trajectories.txt'.                  *
12  !                                                                            *
13  !     Author: A. Stohl                                                       *
14  !                                                                            *
15  !     24 January 2002                                                        *
16  !                                                                            *
17  ! Variables:                                                                 *
18  ! fclust          fraction of particles belonging to each cluster            *
19  ! hmixcenter      mean mixing height for all particles                       *
20  ! ncluster        number of clusters to be used                              *
21  ! pvcenter        mean PV for all particles                                  *
22  ! pvfract         fraction of particles with PV<2pvu                         *
23  ! rms             total horizontal rms distance after clustering             *
24  ! rmsdist         total horizontal rms distance before clustering            *
25  ! rmsclust        horizontal rms distance for each individual cluster        *
26  ! topocenter      mean topography underlying all particles                   *
27  ! tropocenter     mean tropopause height at the positions of particles       *
28  ! tropofract      fraction of particles within the troposphere               *
29  ! zrms            total vertical rms distance after clustering               *
30  ! zrmsdist        total vertical rms distance before clustering              *
31  ! xclust,yclust,  Cluster centroid positions                                 *
32  ! zclust                                                                     *
33  !                                                                            *
34  !*****************************************************************************
35
36  use point_mod
37  use par_mod
38  use com_mod
39  use mean_mod
40
41  implicit none
42
43  integer :: itime,ix,jy,ixp,jyp,indexh,i,j,k,m,n,il,ind,indz,indzp
44  real :: xl(maxpart),yl(maxpart),zl(maxpart)
45  real :: xcenter,ycenter,zcenter,dist,distance,rmsdist,zrmsdist
46
47  real :: xclust(ncluster),yclust(ncluster),zclust(ncluster)
48  real :: fclust(ncluster),rms,rmsclust(ncluster),zrms
49
50  real :: dt1,dt2,dtt,ddx,ddy,rddx,rddy,p1,p2,p3,p4,dz1,dz2,dz
51  real :: topo,topocenter,hm(2),hmixi,hmixfract,hmixcenter
52  real :: pv1(2),pvprof(2),pvi,pvcenter,pvfract,tr(2),tri,tropofract
53  real :: tropocenter
54
55
56  dt1=real(itime-memtime(1))
57  dt2=real(memtime(2)-itime)
58  dtt=1./(dt1+dt2)
59
60
61  ! Loop about all release points
62  !******************************
63
64  do j=1,numpoint
65    if (abs(ireleasestart(j)-itime).gt.lage(nageclass)) goto 10
66    topocenter=0.
67    hmixcenter=0.
68    hmixfract=0.
69    tropocenter=0.
70    tropofract=0.
71    pvfract=0.
72    pvcenter=0.
73    rmsdist=0.
74    zrmsdist=0.
75
76    n=0
77    do i=1,numpart
78      if (itra1(i).ne.itime) goto 20
79      if (npoint(i).ne.j) goto 20
80      n=n+1
81      xl(n)=xlon0+xtra1(i)*dx
82      yl(n)=ylat0+ytra1(i)*dy
83      zl(n)=ztra1(i)
84
85
86  ! Interpolate PBL height, PV, and tropopause height to each
87  ! particle position in order to determine fraction of particles
88  ! within the PBL, above tropopause height, and average PV.
89  ! Interpolate topography, too, and convert to altitude asl
90  !**************************************************************
91
92      ix=int(xtra1(i))
93      jy=int(ytra1(i))
94      ixp=ix+1
95      jyp=jy+1
96      ddx=xtra1(i)-real(ix)
97      ddy=ytra1(i)-real(jy)
98      rddx=1.-ddx
99      rddy=1.-ddy
100      p1=rddx*rddy
101      p2=ddx*rddy
102      p3=rddx*ddy
103      p4=ddx*ddy
104
105  ! Topography
106  !***********
107
108      topo=p1*oro(ix ,jy) &
109           + p2*oro(ixp,jy) &
110           + p3*oro(ix ,jyp) &
111           + p4*oro(ixp,jyp)
112      topocenter=topocenter+topo
113
114  ! Potential vorticity
115  !********************
116
117      do il=2,nz
118        if (height(il).gt.zl(n)) then
119          indz=il-1
120          indzp=il
121          goto 6
122        endif
123      end do
1246     continue
125
126      dz1=zl(n)-height(indz)
127      dz2=height(indzp)-zl(n)
128      dz=1./(dz1+dz2)
129
130
131      do ind=indz,indzp
132        do m=1,2
133          indexh=memind(m)
134          pv1(m)=p1*pv(ix ,jy ,ind,indexh) &
135               +p2*pv(ixp,jy ,ind,indexh) &
136               +p3*pv(ix ,jyp,ind,indexh) &
137               +p4*pv(ixp,jyp,ind,indexh)
138        end do
139        pvprof(ind-indz+1)=(pv1(1)*dt2+pv1(2)*dt1)*dtt
140      end do
141      pvi=(dz1*pvprof(2)+dz2*pvprof(1))*dz
142      pvcenter=pvcenter+pvi
143      if (yl(n).gt.0.) then
144        if (pvi.lt.2.) pvfract=pvfract+1.
145      else
146        if (pvi.gt.-2.) pvfract=pvfract+1.
147      endif
148
149
150  ! Tropopause and PBL height
151  !**************************
152
153      do m=1,2
154        indexh=memind(m)
155
156        tr(m)=p1*tropopause(ix ,jy ,1,indexh) &
157             + p2*tropopause(ixp,jy ,1,indexh) &
158             + p3*tropopause(ix ,jyp,1,indexh) &
159             + p4*tropopause(ixp,jyp,1,indexh)
160
161        hm(m)=p1*hmix(ix ,jy ,1,indexh) &
162             + p2*hmix(ixp,jy ,1,indexh) &
163             + p3*hmix(ix ,jyp,1,indexh) &
164             + p4*hmix(ixp,jyp,1,indexh)
165      end do
166
167      hmixi=(hm(1)*dt2+hm(2)*dt1)*dtt
168      tri=(tr(1)*dt2+tr(2)*dt1)*dtt
169      if (zl(n).lt.tri) tropofract=tropofract+1.
170      tropocenter=tropocenter+tri+topo
171      if (zl(n).lt.hmixi) hmixfract=hmixfract+1.
172      zl(n)=zl(n)+topo        ! convert to height asl
173      hmixcenter=hmixcenter+hmixi
174
175
17620    continue
177    end do
178
179
180  ! Make statistics for all plumes with n>0 particles
181  !**************************************************
182
183    if (n.gt.0) then
184      topocenter=topocenter/real(n)
185      hmixcenter=hmixcenter/real(n)
186      pvcenter=pvcenter/real(n)
187      tropocenter=tropocenter/real(n)
188      hmixfract=100.*hmixfract/real(n)
189      pvfract=100.*pvfract/real(n)
190      tropofract=100.*tropofract/real(n)
191
192  ! Cluster the particle positions
193  !*******************************
194
195      call clustering(xl,yl,zl,n,xclust,yclust,zclust,fclust,rms, &
196           rmsclust,zrms)
197
198
199  ! Determine center of mass position on earth and average height
200  !**************************************************************
201
202      call centerofmass(xl,yl,n,xcenter,ycenter)
203      call mean(zl,zcenter,zrmsdist,n)
204
205  ! Root mean square distance from center of mass
206  !**********************************************
207
208      do k=1,n
209        dist=distance(yl(k),xl(k),ycenter,xcenter)
210        rmsdist=rmsdist+dist*dist
211      end do
212      if (rmsdist.gt.0.) rmsdist=sqrt(rmsdist/real(n))
213      rmsdist=max(rmsdist,0.)
214
215  ! Write out results in trajectory data file
216  !******************************************
217
218      write(unitouttraj,'(i5,i8,2f9.4,4f8.1,f8.2,4f8.1,3f6.1,&
219           &5(2f8.3,f7.0,f6.1,f8.1))')&
220           &j,itime-(ireleasestart(j)+ireleaseend(j))/2, &
221           xcenter,ycenter,zcenter,topocenter,hmixcenter,tropocenter, &
222           pvcenter,rmsdist,rms,zrmsdist,zrms,hmixfract,pvfract, &
223           tropofract, &
224           (xclust(k),yclust(k),zclust(k),fclust(k),rmsclust(k), &
225           k=1,ncluster)
226    endif
227
228
22910  continue
230  end do
231
232
233end subroutine plumetraj
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG