[e200b7a] | 1 | subroutine clustering(xl,yl,zl,n,xclust,yclust,zclust,fclust,rms, & |
---|
| 2 | rmsclust,zrms) |
---|
| 3 | ! i i i i o o o o o |
---|
| 4 | ! o o |
---|
| 5 | !***************************************************************************** |
---|
| 6 | ! * |
---|
| 7 | ! This routine clusters the particle position into ncluster custers. * |
---|
| 8 | ! Input are the longitudes (xl) and latitudes (yl) of the individual * |
---|
| 9 | ! points, output are the cluster mean positions (xclust,yclust). * |
---|
| 10 | ! Vertical positions are not directly used for the clustering. * |
---|
| 11 | ! * |
---|
| 12 | ! For clustering, the procedure described in Dorling et al. (1992) is used.* |
---|
| 13 | ! * |
---|
| 14 | ! Dorling, S.R., Davies, T.D. and Pierce, C.E. (1992): * |
---|
| 15 | ! Cluster analysis: a technique for estimating the synoptic meteorological * |
---|
| 16 | ! controls on air and precipitation chemistry - method and applications. * |
---|
| 17 | ! Atmospheric Environment 26A, 2575-2581. * |
---|
| 18 | ! * |
---|
| 19 | ! * |
---|
| 20 | ! Author: A. Stohl * |
---|
| 21 | ! * |
---|
| 22 | ! 1 February 2002 * |
---|
| 23 | ! * |
---|
| 24 | ! Variables: * |
---|
| 25 | ! fclust fraction of particles belonging to each cluster * |
---|
| 26 | ! ncluster number of clusters to be used * |
---|
| 27 | ! rms total horizontal rms distance after clustering * |
---|
| 28 | ! rmsclust horizontal rms distance for each individual cluster * |
---|
| 29 | ! zrms total vertical rms distance after clustering * |
---|
| 30 | ! xclust,yclust, Cluster centroid positions * |
---|
| 31 | ! zclust * |
---|
| 32 | ! xl,yl,zl particle positions * |
---|
| 33 | ! * |
---|
| 34 | !***************************************************************************** |
---|
| 35 | |
---|
| 36 | use par_mod |
---|
| 37 | |
---|
| 38 | implicit none |
---|
| 39 | |
---|
| 40 | integer :: n,i,j,l,nclust(maxpart),numb(ncluster),ncl |
---|
| 41 | real :: xl(n),yl(n),zl(n),xclust(ncluster),yclust(ncluster),x,y,z |
---|
| 42 | real :: zclust(ncluster),distance2,distances,distancemin,rms,rmsold |
---|
| 43 | real :: xav(ncluster),yav(ncluster),zav(ncluster),fclust(ncluster) |
---|
| 44 | real :: rmsclust(ncluster) |
---|
| 45 | real :: zdist,zrms |
---|
| 46 | |
---|
| 47 | |
---|
| 48 | |
---|
| 49 | if (n.lt.ncluster) return |
---|
| 50 | rmsold=-5. |
---|
| 51 | |
---|
| 52 | ! Convert longitude and latitude from degrees to radians |
---|
| 53 | !******************************************************* |
---|
| 54 | |
---|
| 55 | do i=1,n |
---|
| 56 | nclust(i)=i |
---|
| 57 | xl(i)=xl(i)*pi180 |
---|
| 58 | yl(i)=yl(i)*pi180 |
---|
| 59 | end do |
---|
| 60 | |
---|
| 61 | |
---|
| 62 | ! Generate a seed for each cluster |
---|
| 63 | !********************************* |
---|
| 64 | |
---|
| 65 | do j=1,ncluster |
---|
| 66 | zclust(j)=0. |
---|
| 67 | xclust(j)=xl(j*n/ncluster) |
---|
| 68 | yclust(j)=yl(j*n/ncluster) |
---|
| 69 | end do |
---|
| 70 | |
---|
| 71 | |
---|
| 72 | ! Iterative loop to compute the cluster means |
---|
| 73 | !******************************************** |
---|
| 74 | |
---|
| 75 | do l=1,100 |
---|
| 76 | |
---|
| 77 | ! Assign each particle to a cluster: criterion minimum distance to the |
---|
| 78 | ! cluster mean position |
---|
| 79 | !********************************************************************* |
---|
| 80 | |
---|
| 81 | |
---|
| 82 | do i=1,n |
---|
| 83 | distancemin=10.**10. |
---|
| 84 | do j=1,ncluster |
---|
| 85 | distances=distance2(yl(i),xl(i),yclust(j),xclust(j)) |
---|
| 86 | if (distances.lt.distancemin) then |
---|
| 87 | distancemin=distances |
---|
| 88 | ncl=j |
---|
| 89 | endif |
---|
| 90 | end do |
---|
| 91 | nclust(i)=ncl |
---|
| 92 | end do |
---|
| 93 | |
---|
| 94 | |
---|
| 95 | ! Recalculate the cluster centroid position: convert to 3D Cartesian coordinates, |
---|
| 96 | ! calculate mean position, and re-project this point onto the Earth's surface |
---|
| 97 | !***************************************************************************** |
---|
| 98 | |
---|
| 99 | do j=1,ncluster |
---|
| 100 | xav(j)=0. |
---|
| 101 | yav(j)=0. |
---|
| 102 | zav(j)=0. |
---|
| 103 | rmsclust(j)=0. |
---|
| 104 | numb(j)=0 |
---|
| 105 | end do |
---|
| 106 | rms=0. |
---|
| 107 | |
---|
| 108 | do i=1,n |
---|
| 109 | numb(nclust(i))=numb(nclust(i))+1 |
---|
| 110 | distances=distance2(yl(i),xl(i), & |
---|
| 111 | yclust(nclust(i)),xclust(nclust(i))) |
---|
| 112 | |
---|
| 113 | ! rms is the total rms of all particles |
---|
| 114 | ! rmsclust is the rms for a particular cluster |
---|
| 115 | !********************************************* |
---|
| 116 | |
---|
| 117 | rms=rms+distances*distances |
---|
| 118 | rmsclust(nclust(i))=rmsclust(nclust(i))+distances*distances |
---|
| 119 | |
---|
| 120 | ! Calculate Cartesian 3D coordinates from longitude and latitude |
---|
| 121 | !*************************************************************** |
---|
| 122 | |
---|
| 123 | x = cos(yl(i))*sin(xl(i)) |
---|
| 124 | y = -1.*cos(yl(i))*cos(xl(i)) |
---|
| 125 | z = sin(yl(i)) |
---|
| 126 | xav(nclust(i))=xav(nclust(i))+x |
---|
| 127 | yav(nclust(i))=yav(nclust(i))+y |
---|
| 128 | zav(nclust(i))=zav(nclust(i))+z |
---|
| 129 | end do |
---|
| 130 | |
---|
| 131 | rms=sqrt(rms/real(n)) |
---|
| 132 | |
---|
| 133 | |
---|
| 134 | ! Find the mean location in Cartesian coordinates |
---|
| 135 | !************************************************ |
---|
| 136 | |
---|
| 137 | do j=1,ncluster |
---|
| 138 | if (numb(j).gt.0) then |
---|
| 139 | rmsclust(j)=sqrt(rmsclust(j)/real(numb(j))) |
---|
| 140 | xav(j)=xav(j)/real(numb(j)) |
---|
| 141 | yav(j)=yav(j)/real(numb(j)) |
---|
| 142 | zav(j)=zav(j)/real(numb(j)) |
---|
| 143 | |
---|
| 144 | ! Project the point back onto Earth's surface |
---|
| 145 | !******************************************** |
---|
| 146 | |
---|
| 147 | xclust(j)=atan2(xav(j),-1.*yav(j)) |
---|
| 148 | yclust(j)=atan2(zav(j),sqrt(xav(j)*xav(j)+yav(j)*yav(j))) |
---|
| 149 | endif |
---|
| 150 | end do |
---|
| 151 | |
---|
| 152 | |
---|
| 153 | ! Leave the loop if the RMS distance decreases only slightly between 2 iterations |
---|
| 154 | !***************************************************************************** |
---|
| 155 | |
---|
| 156 | if ((l.gt.1).and.(abs(rms-rmsold)/rmsold.lt.0.005)) goto 99 |
---|
| 157 | rmsold=rms |
---|
| 158 | |
---|
| 159 | end do |
---|
| 160 | |
---|
| 161 | 99 continue |
---|
| 162 | |
---|
| 163 | ! Convert longitude and latitude from radians to degrees |
---|
| 164 | !******************************************************* |
---|
| 165 | |
---|
| 166 | do i=1,n |
---|
| 167 | xl(i)=xl(i)/pi180 |
---|
| 168 | yl(i)=yl(i)/pi180 |
---|
| 169 | zclust(nclust(i))=zclust(nclust(i))+zl(i) |
---|
| 170 | end do |
---|
| 171 | |
---|
| 172 | do j=1,ncluster |
---|
| 173 | xclust(j)=xclust(j)/pi180 |
---|
| 174 | yclust(j)=yclust(j)/pi180 |
---|
| 175 | if (numb(j).gt.0) zclust(j)=zclust(j)/real(numb(j)) |
---|
| 176 | fclust(j)=100.*real(numb(j))/real(n) |
---|
| 177 | end do |
---|
| 178 | |
---|
| 179 | ! Determine total vertical RMS deviation |
---|
| 180 | !*************************************** |
---|
| 181 | |
---|
| 182 | zrms=0. |
---|
| 183 | do i=1,n |
---|
| 184 | zdist=zl(i)-zclust(nclust(i)) |
---|
| 185 | zrms=zrms+zdist*zdist |
---|
| 186 | end do |
---|
| 187 | if (zrms.gt.0.) zrms=sqrt(zrms/real(n)) |
---|
| 188 | |
---|
| 189 | end subroutine clustering |
---|