CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/RecoLocalMuon/CSCSegment/src/CSCSegAlgoPreClustering.cc

Go to the documentation of this file.
00001 
00013 #include "RecoLocalMuon/CSCSegment/src/CSCSegAlgoPreClustering.h"
00014 
00015 #include "Geometry/CSCGeometry/interface/CSCLayer.h"
00016 #include <Geometry/CSCGeometry/interface/CSCChamber.h>
00017 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00018 
00019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00021 
00022 #include <algorithm>
00023 #include <cmath>
00024 #include <iostream>
00025 #include <string>
00026 
00027 
00028 /* Constructor
00029  *
00030  */
00031 CSCSegAlgoPreClustering::CSCSegAlgoPreClustering(const edm::ParameterSet& ps) {
00032   dXclusBoxMax           = ps.getParameter<double>("dXclusBoxMax");
00033   dYclusBoxMax           = ps.getParameter<double>("dYclusBoxMax");
00034   debug                  = ps.getUntrackedParameter<bool>("CSCSegmentDebug");
00035 }
00036 
00037 
00038 /* Destructor:
00039  *
00040  */
00041 CSCSegAlgoPreClustering::~CSCSegAlgoPreClustering(){
00042 
00043 }
00044 
00045 
00046 /* clusterHits
00047  *
00048  */
00049 std::vector< std::vector<const CSCRecHit2D*> > 
00050 CSCSegAlgoPreClustering::clusterHits( const CSCChamber* aChamber, ChamberHitContainer rechits) {
00051 
00052   theChamber = aChamber;
00053 
00054   std::vector<ChamberHitContainer> rechits_clusters; // this is a collection of groups of rechits
00055 
00056   //float dXclus = 0.0;
00057   //float dYclus = 0.0;
00058   float dXclus_box = 0.0;
00059   float dYclus_box = 0.0;
00060 
00061   std::vector<const CSCRecHit2D*> temp;
00062 
00063   std::vector< ChamberHitContainer > seeds;
00064 
00065   std::vector<float> running_meanX;
00066   std::vector<float> running_meanY;
00067 
00068   std::vector<float> seed_minX;
00069   std::vector<float> seed_maxX;
00070   std::vector<float> seed_minY;
00071   std::vector<float> seed_maxY;
00072 
00073   // split rechits into subvectors and return vector of vectors:
00074   // Loop over rechits 
00075   // Create one seed per hit
00076     for(unsigned int i = 0; i < rechits.size(); ++i) {
00077 
00078         temp.clear();
00079 
00080         temp.push_back(rechits[i]);
00081 
00082         seeds.push_back(temp);
00083 
00084         // First added hit in seed defines the mean to which the next hit is compared
00085         // for this seed.
00086 
00087         running_meanX.push_back( rechits[i]->localPosition().x() );
00088         running_meanY.push_back( rechits[i]->localPosition().y() );
00089         
00090         // set min/max X and Y for box containing the hits in the precluster:
00091         seed_minX.push_back( rechits[i]->localPosition().x() );
00092         seed_maxX.push_back( rechits[i]->localPosition().x() );
00093         seed_minY.push_back( rechits[i]->localPosition().y() );
00094         seed_maxY.push_back( rechits[i]->localPosition().y() );
00095     }
00096     
00097     // merge clusters that are too close
00098     // measure distance between final "running mean"
00099       for(size_t NNN = 0; NNN < seeds.size(); ++NNN) {
00100         
00101         for(size_t MMM = NNN+1; MMM < seeds.size(); ++MMM) {
00102           if(running_meanX[MMM] == 999999. || running_meanX[NNN] == 999999. ) {
00103             std::cout<<"We should never see this line now!!!"<<std::endl;
00104             continue; //skip seeds that have been used 
00105           }
00106           
00107           // calculate cut criteria for simple running mean distance cut:
00108           //dXclus = fabs(running_meanX[NNN] - running_meanX[MMM]);
00109           //dYclus = fabs(running_meanY[NNN] - running_meanY[MMM]);
00110 
00111           // calculate minmal distance between precluster boxes containing the hits:
00112           if ( running_meanX[NNN] > running_meanX[MMM] ) dXclus_box = seed_minX[NNN] - seed_maxX[MMM];
00113           else                                           dXclus_box = seed_minX[MMM] - seed_maxX[NNN];
00114           if ( running_meanY[NNN] > running_meanY[MMM] ) dYclus_box = seed_minY[NNN] - seed_maxY[MMM];
00115           else                                           dYclus_box = seed_minY[MMM] - seed_maxY[NNN];
00116           
00117           
00118           if( dXclus_box < dXclusBoxMax && dYclus_box < dYclusBoxMax ) {
00119             // merge clusters!
00120             // merge by adding seed NNN to seed MMM and erasing seed NNN
00121             
00122             // calculate running mean for the merged seed:
00123             running_meanX[MMM] = (running_meanX[NNN]*seeds[NNN].size() + running_meanX[MMM]*seeds[MMM].size()) / (seeds[NNN].size()+seeds[MMM].size());
00124             running_meanY[MMM] = (running_meanY[NNN]*seeds[NNN].size() + running_meanY[MMM]*seeds[MMM].size()) / (seeds[NNN].size()+seeds[MMM].size());
00125             
00126             // update min/max X and Y for box containing the hits in the merged cluster:
00127             if ( seed_minX[NNN] <= seed_minX[MMM] ) seed_minX[MMM] = seed_minX[NNN];
00128             if ( seed_maxX[NNN] >  seed_maxX[MMM] ) seed_maxX[MMM] = seed_maxX[NNN];
00129             if ( seed_minY[NNN] <= seed_minY[MMM] ) seed_minY[MMM] = seed_minY[NNN];
00130             if ( seed_maxY[NNN] >  seed_maxY[MMM] ) seed_maxY[MMM] = seed_maxY[NNN];
00131             
00132             // add seed NNN to MMM (lower to larger number)
00133             seeds[MMM].insert(seeds[MMM].end(),seeds[NNN].begin(),seeds[NNN].end());
00134             
00135             // mark seed NNN as used (at the moment just set running mean to 999999.)
00136             running_meanX[NNN] = 999999.;
00137             running_meanY[NNN] = 999999.;
00138             // we have merged a seed (NNN) to the highter seed (MMM) - need to contimue to 
00139             // next seed (NNN+1)
00140             break;
00141           }
00142 
00143         }
00144       }
00145 
00146       // hand over the final seeds to the output
00147       // would be more elegant if we could do the above step with 
00148       // erasing the merged ones, rather than the 
00149       for(size_t NNN = 0; NNN < seeds.size(); ++NNN) {
00150         if (running_meanX[NNN] == 999999.) continue; //skip seeds that have been marked as used up in merging
00151         rechits_clusters.push_back(seeds[NNN]);
00152         mean_x = running_meanX[NNN];
00153         mean_y = running_meanY[NNN];
00154         err_x  = (seed_maxX[NNN]-seed_minX[NNN])/3.464101615; // use box size divided by sqrt(12) as position error estimate
00155         err_y  = (seed_maxY[NNN]-seed_minY[NNN])/3.464101615; // use box size divided by sqrt(12) as position error estimate
00156 
00157       }
00158 
00159   return rechits_clusters; 
00160 }