CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/RecoMuon/MuonSeedGenerator/src/SETPatternRecognition.cc

Go to the documentation of this file.
00001 
00005 #include "RecoMuon/MuonSeedGenerator/src/SETPatternRecognition.h"
00006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00007 #include "FWCore/Framework/interface/Event.h"
00008 #include "RecoMuon/TrackingTools/interface/MuonPatternRecoDumper.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00011 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00012 #include "DataFormats/DTRecHit/interface/DTRecSegment4D.h"
00013 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"            
00014 #include "DataFormats/CSCRecHit/interface/CSCSegmentCollection.h"             
00015 #include "DataFormats/RPCRecHit/interface/RPCRecHitCollection.h"                
00016 #include "TMath.h"
00017 
00018 
00019 using namespace edm;
00020 using namespace std;
00021 
00022 SETPatternRecognition::SETPatternRecognition(const ParameterSet& parameterSet)
00023 : MuonSeedVPatternRecognition(parameterSet.getParameter<ParameterSet>("SETTrajBuilderParameters").getParameter<ParameterSet>("FilterParameters"))
00024 {
00025   const string metname = "Muon|RecoMuon|SETPatternRecognition";  
00026   // Parameter set for the Builder
00027   ParameterSet trajectoryBuilderParameters = parameterSet.getParameter<ParameterSet>("SETTrajBuilderParameters");
00028   // The inward-outward fitter (starts from seed state)
00029   ParameterSet filterPSet = trajectoryBuilderParameters.getParameter<ParameterSet>("FilterParameters");
00030   maxActiveChambers = filterPSet.getParameter<int>("maxActiveChambers"); 
00031   useRPCs = filterPSet.getParameter<bool>("EnableRPCMeasurement");
00032   DTRecSegmentLabel  = filterPSet.getParameter<edm::InputTag>("DTRecSegmentLabel");
00033   CSCRecSegmentLabel  = filterPSet.getParameter<edm::InputTag>("CSCRecSegmentLabel");
00034   RPCRecSegmentLabel  = filterPSet.getParameter<edm::InputTag>("RPCRecSegmentLabel");
00035 
00036   outsideChamberErrorScale = filterPSet.getParameter<double>("OutsideChamberErrorScale");
00037   minLocalSegmentAngle = filterPSet.getParameter<double>("MinLocalSegmentAngle");
00038   //----
00039 
00040 } 
00041 
00042 //---- it is a "cluster recognition" actually; the pattern recognition is in SETFilter 
00043 void SETPatternRecognition::produce(const edm::Event& event, const edm::EventSetup& eventSetup,
00044                          std::vector<MuonRecHitContainer> & segments_clusters)
00045 {
00046   const string metname = "Muon|RecoMuon|SETMuonSeedSeed";  
00047 
00048   //---- Build collection of all segments
00049   MuonRecHitContainer muonRecHits;
00050   MuonRecHitContainer muonRecHits_DT2D_hasPhi;
00051   MuonRecHitContainer muonRecHits_DT2D_hasZed;
00052   MuonRecHitContainer muonRecHits_RPC;
00053 
00054   // ********************************************;
00055   // Get the DT-Segment collection from the Event
00056   // ********************************************;
00057 
00058   edm::Handle<DTRecSegment4DCollection> dtRecHits;
00059   event.getByLabel(DTRecSegmentLabel, dtRecHits);
00060   std::vector<DTChamberId> chambers_DT;
00061   std::vector<DTChamberId>::const_iterator chIt_DT;
00062   for (DTRecSegment4DCollection::const_iterator rechit = dtRecHits->begin(); rechit!=dtRecHits->end();++rechit) {
00063     bool insert = true;
00064     for(chIt_DT=chambers_DT.begin(); chIt_DT != chambers_DT.end(); ++chIt_DT){
00065       if (
00066           ((*rechit).chamberId().wheel()) == ((*chIt_DT).wheel()) &&
00067           ((*rechit).chamberId().station() == (*chIt_DT).station()) &&
00068           ((*rechit).chamberId().sector() == (*chIt_DT).sector())){
00069         insert = false;
00070       }
00071     }
00072     if (insert){
00073       chambers_DT.push_back((*rechit).chamberId());
00074     }
00075     if(segmentCleaning((*rechit).geographicalId(), 
00076                        rechit->localPosition(), rechit->localPositionError(),
00077                        rechit->localDirection(), rechit->localDirectionError(),
00078                        rechit->chi2(), rechit->degreesOfFreedom())){
00079       continue;
00080     }
00081     if( (rechit->hasZed() && rechit->hasPhi()) ) {
00082     muonRecHits.push_back(MuonTransientTrackingRecHit::specificBuild(theService->trackingGeometry()->idToDet((*rechit).geographicalId()),&*rechit));
00083     }
00084     else if(rechit->hasZed()) {
00085     muonRecHits_DT2D_hasZed.push_back(MuonTransientTrackingRecHit::specificBuild(theService->trackingGeometry()->idToDet((*rechit).geographicalId()),&*rechit));
00086     }
00087     else if(rechit->hasPhi()) { // safeguard
00088     muonRecHits_DT2D_hasPhi.push_back(MuonTransientTrackingRecHit::specificBuild(theService->trackingGeometry()->idToDet((*rechit).geographicalId()),&*rechit));
00089     }
00090     else {
00091       //std::cout<<"Warning in "<<metname<<": DT segment which claims to have neither phi nor Z."<<std::endl;
00092     }
00093   }
00094   //std::cout<<"DT done"<<std::endl;
00095 
00096   // ********************************************;
00097   // Get the CSC-Segment collection from the event
00098   // ********************************************;
00099 
00100   edm::Handle<CSCSegmentCollection> cscSegments;
00101   event.getByLabel(CSCRecSegmentLabel, cscSegments);
00102   std::vector<CSCDetId> chambers_CSC;
00103   std::vector<CSCDetId>::const_iterator chIt_CSC;
00104   for(CSCSegmentCollection::const_iterator rechit=cscSegments->begin(); rechit != cscSegments->end(); ++rechit) {
00105     bool insert = true;
00106     for(chIt_CSC=chambers_CSC.begin(); chIt_CSC != chambers_CSC.end(); ++chIt_CSC){
00107       if (((*rechit).cscDetId().chamber() == (*chIt_CSC).chamber()) &&
00108           ((*rechit).cscDetId().station() == (*chIt_CSC).station()) &&
00109           ((*rechit).cscDetId().ring() == (*chIt_CSC).ring()) &&
00110           ((*rechit).cscDetId().endcap() == (*chIt_CSC).endcap())){
00111         insert = false;
00112       }
00113     }
00114     if (insert){
00115       chambers_CSC.push_back((*rechit).cscDetId().chamberId());
00116     }
00117     if(segmentCleaning((*rechit).geographicalId(), 
00118                        rechit->localPosition(), rechit->localPositionError(),
00119                        rechit->localDirection(), rechit->localDirectionError(),
00120                        rechit->chi2(), rechit->degreesOfFreedom())){
00121       continue;
00122     }
00123     muonRecHits.push_back(MuonTransientTrackingRecHit::specificBuild(theService->trackingGeometry()->idToDet((*rechit).geographicalId()),&*rechit));
00124   }
00125   //std::cout<<"CSC done"<<std::endl;
00126 
00127   // ********************************************;
00128   // Get the RPC-Hit collection from the event
00129   // ********************************************;
00130 
00131   edm::Handle<RPCRecHitCollection> rpcRecHits;
00132   event.getByLabel(RPCRecSegmentLabel, rpcRecHits);
00133   if(useRPCs){
00134     for(RPCRecHitCollection::const_iterator rechit=rpcRecHits->begin(); rechit != rpcRecHits->end(); ++rechit) {
00135       // RPCs are special
00136       const LocalVector  localDirection(0.,0.,1.);
00137       const LocalError localDirectionError (0.,0.,0.); 
00138       const double chi2 = 1.;
00139       const int ndf = 1;
00140       if(segmentCleaning((*rechit).geographicalId(), 
00141                          rechit->localPosition(), rechit->localPositionError(),
00142                          localDirection, localDirectionError,
00143                          chi2, ndf)){
00144         continue;
00145       }
00146       muonRecHits_RPC.push_back(MuonTransientTrackingRecHit::specificBuild(theService->trackingGeometry()->idToDet((*rechit).geographicalId()),&*rechit));
00147     }
00148   }
00149   //std::cout<<"RPC done"<<std::endl;
00150   //
00151   if(int(chambers_DT.size() + chambers_CSC.size()) > maxActiveChambers){
00152     // std::cout <<" Too many active chambers : nDT = "<<chambers_DT.size()<<
00153     // " nCSC = "<<chambers_CSC.size()<<"  Skip them all."<<std::endl;
00154     edm::LogWarning("tooManyActiveChambers")<<" Too many active chambers : nDT = "<<chambers_DT.size()
00155                      <<" nCSC = "<<chambers_CSC.size()<<"  Skip them all.";
00156     muonRecHits.clear();                                
00157     muonRecHits_DT2D_hasPhi.clear();    
00158     muonRecHits_DT2D_hasZed.clear();
00159     muonRecHits_RPC.clear();
00160   }
00161   //---- Find "pre-clusters" from all segments; these contain potential muon candidates
00162 
00163   //---- From all the hits (i.e. segments; sometimes "rechits" is also used with the same meaning;
00164   //---- this convention has meaning in the global reconstruction though could be misleading 
00165   //---- from a local reconstruction point of view; "local rechits" are used in the backward fit only) 
00166   //---- make clusters of hits; a trajectory could contain hits from one cluster only   
00167 
00168   // the clustering procedure is very similar to the one used in the segment reconstruction 
00169 
00170   bool useDT2D_hasPhi = true;
00171   bool useDT2D_hasZed = true;
00172   double dXclusBoxMax         = 0.60; // phi - can be as large as 15 - 20 degrees for 6 GeV muons
00173   double dYclusBoxMax = 0.;
00174 
00175   // this is the main selection criteria; the value of 0.02 rad seems wide enough to 
00176   // contain any hit from a passing muon and still narrow enough to remove good part of
00177   // possible "junk" hits
00178   // (Comment: it may be better to allow maximum difference between any two hits in a trajectory
00179   // to be 0.02 or 0.04 or ...; currently the requirement below is imposed on two consecutive hits)  
00180 
00181   dYclusBoxMax              = 0.02; // theta // hardoded - remove it!
00182   
00183   // X and Y are distance variables - we use eta and phi here
00184 
00185   float dXclus = 0.0;
00186   float dXclus_box = 0.0;
00187   float dYclus_box = 0.0;
00188 
00189   MuonRecHitContainer temp;
00190 
00191   std::vector< MuonRecHitContainer > seeds;
00192 
00193   std::vector<float> running_meanX;
00194   std::vector<float> running_meanY;
00195 
00196   std::vector<float> seed_minX;
00197   std::vector<float> seed_maxX;
00198   std::vector<float> seed_minY;
00199   std::vector<float> seed_maxY;
00200 
00201   // split rechits into subvectors and return vector of vectors:
00202   // Loop over rechits 
00203   // Create one seed per hit
00204   for (MuonRecHitContainer::const_iterator it = muonRecHits.begin(); it != muonRecHits.end(); ++it ) {
00205 
00206     // try to avoid using 2D DT segments. We will add them later to the 
00207     // clusters they are most likely to belong to. Might need to add them 
00208     // to more than just one cluster, if we find them to be consistent with 
00209     // more than one. This would lead to an implicit sharing of hits between 
00210     // SA muon candidates. 
00211 
00212     temp.clear();
00213 
00214     temp.push_back((*it));
00215 
00216     seeds.push_back(temp);
00217 
00218     // First added hit in seed defines the mean to which the next hit is compared
00219     // for this seed.
00220 
00221     running_meanX.push_back( (*it)->globalPosition().phi() );
00222     running_meanY.push_back( (*it)->globalPosition().theta() );
00223 
00224     // set min/max X and Y for box containing the hits in the precluster:
00225     seed_minX.push_back( (*it)->globalPosition().phi() );
00226     seed_maxX.push_back( (*it)->globalPosition().phi() );
00227     seed_minY.push_back( (*it)->globalPosition().theta() );
00228     seed_maxY.push_back( (*it)->globalPosition().theta() );
00229   }
00230 
00231   // merge clusters that are too close
00232   // measure distance between final "running mean"
00233   for(unsigned int NNN = 0; NNN < seeds.size(); ++NNN) {
00234 
00235     for(unsigned int MMM = NNN+1; MMM < seeds.size(); ++MMM) {
00236       if(running_meanX[MMM] == 999999. || running_meanX[NNN] == 999999. ) {
00237         //        LogDebug("CSC") << "CSCSegmentST::clusterHits: Warning: Skipping used seeds, this should happen - inform developers!\n";
00238         //std::cout<<"We should never see this line now!!!"<<std::endl;
00239         continue; //skip seeds that have been used 
00240       }
00241 
00242       // Some complications for using phi as a clustering variable due to wrap-around (-pi = pi)
00243       // Define temporary mean, min, and max variables for the cluster which could be merged (NNN)
00244       double temp_meanX = running_meanX[NNN];
00245       double temp_minX = seed_minX[NNN];
00246       double temp_maxX = seed_maxX[NNN];
00247 
00248       // check if the difference between the two phi values is greater than pi
00249       // if so, need to shift temporary values by 2*pi to make a valid comparison
00250       dXclus = running_meanX[NNN] - running_meanX[MMM];
00251       if (dXclus > TMath::Pi()) {
00252         temp_meanX = temp_meanX - 2.*TMath::Pi();
00253         temp_minX = temp_minX - 2.*TMath::Pi();
00254         temp_maxX = temp_maxX - 2.*TMath::Pi();
00255       }
00256       if (dXclus < -TMath::Pi()) {
00257         temp_meanX = temp_meanX + 2.*TMath::Pi();
00258         temp_minX = temp_minX + 2.*TMath::Pi();
00259         temp_maxX = temp_maxX + 2.*TMath::Pi();
00260       }
00261 
00262       //       // calculate cut criteria for simple running mean distance cut:
00263       //       // not sure that these values are really used anywhere
00264 
00265       // calculate minmal distance between precluster boxes containing the hits:
00266       // use the temp variables from above for phi of the NNN cluster 
00267       if ( temp_meanX > running_meanX[MMM] )         dXclus_box = temp_minX - seed_maxX[MMM];
00268       else                                           dXclus_box = seed_minX[MMM] - temp_maxX;
00269       if ( running_meanY[NNN] > running_meanY[MMM] ) dYclus_box = seed_minY[NNN] - seed_maxY[MMM];
00270       else                                           dYclus_box = seed_minY[MMM] - seed_maxY[NNN];
00271 
00272 
00273       if( dXclus_box < dXclusBoxMax && dYclus_box < dYclusBoxMax ) {
00274         // merge clusters!
00275         // merge by adding seed NNN to seed MMM and erasing seed NNN
00276 
00277         // calculate running mean for the merged seed:
00278         // use the temp variables from above for phi of the NNN cluster 
00279         running_meanX[MMM] = (temp_meanX*seeds[NNN].size() + running_meanX[MMM]*seeds[MMM].size()) / (seeds[NNN].size()+seeds[MMM].size());
00280         running_meanY[MMM] = (running_meanY[NNN]*seeds[NNN].size() + running_meanY[MMM]*seeds[MMM].size()) / (seeds[NNN].size()+seeds[MMM].size());
00281 
00282         // update min/max X and Y for box containing the hits in the merged cluster:
00283         // use the temp variables from above for phi of the NNN cluster 
00284         if ( temp_minX <= seed_minX[MMM] ) seed_minX[MMM] = temp_minX;
00285         if ( temp_maxX >  seed_maxX[MMM] ) seed_maxX[MMM] = temp_maxX;
00286         if ( seed_minY[NNN] <= seed_minY[MMM] ) seed_minY[MMM] = seed_minY[NNN];
00287         if ( seed_maxY[NNN] >  seed_maxY[MMM] ) seed_maxY[MMM] = seed_maxY[NNN];
00288 
00289         // now check to see if the running mean has moved outside of the allowed -pi to pi region
00290         // if so, then adjust shift all values up or down by 2 * pi
00291         if (running_meanX[MMM] > TMath::Pi()) {
00292           running_meanX[MMM] = running_meanX[MMM] - 2.*TMath::Pi();
00293           seed_minX[MMM] = seed_minX[MMM] - 2.*TMath::Pi();
00294           seed_maxX[MMM] = seed_maxX[MMM] - 2.*TMath::Pi();
00295         }
00296         if (running_meanX[MMM] < -TMath::Pi()) {
00297           running_meanX[MMM] = running_meanX[MMM] + 2.*TMath::Pi();
00298           seed_minX[MMM] = seed_minX[MMM] + 2.*TMath::Pi();
00299           seed_maxX[MMM] = seed_maxX[MMM] + 2.*TMath::Pi();
00300         }
00301 
00302         // add seed NNN to MMM (lower to larger number)
00303         seeds[MMM].insert(seeds[MMM].end(),seeds[NNN].begin(),seeds[NNN].end());
00304 
00305         // mark seed NNN as used (at the moment just set running mean to 999999.)
00306         running_meanX[NNN] = 999999.;
00307         running_meanY[NNN] = 999999.;
00308         // we have merged a seed (NNN) to the highter seed (MMM) - need to contimue to 
00309         // next seed (NNN+1)
00310         break;
00311       }
00312 
00313     }
00314   }
00315   bool tooCloseClusters = false;
00316   if(seeds.size()>1){
00317     std::vector <double> seedTheta(seeds.size());
00318     for(unsigned int iSeed = 0;iSeed<seeds.size();++iSeed){
00319       seedTheta[iSeed] = seeds[iSeed][0]->globalPosition().theta(); 
00320       if(iSeed){
00321         double dTheta = fabs(seedTheta[iSeed] - seedTheta[iSeed-1]);
00322         if (dTheta < 0.5){ //? should be something more clever
00323           tooCloseClusters = true;
00324           break;
00325         }
00326       }
00327     }
00328     
00329   }
00330 
00331   // have formed clusters from all hits except for 2D DT segments. Now add the 2D segments to the 
00332   // compatible clusters. For this we compare the mean cluster postition with the 
00333   // 2D segment position. We should use the valid coordinate only and use the bad coordinate 
00334   // as a cross check.
00335   for(unsigned int NNN = 0; NNN < seeds.size(); ++NNN) {
00336     if(running_meanX[NNN] == 999999.) continue; //skip seeds that have been marked as used up in merging
00337 
00338     // We have a valid cluster - loop over all 2D segments.
00339     if(useDT2D_hasZed) {
00340       for (MuonRecHitContainer::const_iterator it2 = muonRecHits_DT2D_hasZed.begin(); it2 != muonRecHits_DT2D_hasZed.end(); ++it2 ) {
00341         // check that global theta of 2-D segment lies within cluster box plus or minus allowed slop
00342         if (((*it2)->globalPosition().theta() < seed_maxY[NNN] + dYclusBoxMax) && ((*it2)->globalPosition().theta() > seed_minY[NNN] - dYclusBoxMax)) {
00343           // check that global phi of 2-D segment (assumed to be center of chamber since no phi hit info)
00344           // matches with cluster box plus or minus allowed slop given that the true phi value could be 
00345           // anywhere within a given chamber (+/- 5 degrees ~ 0.09 radians from center) 
00346           if(
00347              !( 
00348                (
00349                 ((*it2)->globalPosition().phi() + 0.09) < (seed_minX[NNN] - dXclusBoxMax) 
00350                 && 
00351                 ((*it2)->globalPosition().phi() - 0.09) < (seed_minX[NNN] - dXclusBoxMax)
00352                 )
00353                ||
00354                (
00355                 ((*it2)->globalPosition().phi() + 0.09) > (seed_maxX[NNN] + dXclusBoxMax) 
00356                 && 
00357                 ((*it2)->globalPosition().phi() - 0.09) > (seed_maxX[NNN] + dXclusBoxMax)
00358                 )
00359                )
00360              ) { // we have checked that the 2Dsegment is within tight theta boundaries and loose phi boundaries of the current cluster -> add it
00361             seeds[NNN].push_back((*it2));
00362             
00363           }
00364         }                 
00365       }
00366       
00367     }
00368     
00369     // put DT hasphi loop here
00370     if (useDT2D_hasPhi) {
00371       
00372       for (MuonRecHitContainer::const_iterator it2 = muonRecHits_DT2D_hasPhi.begin(); it2 != muonRecHits_DT2D_hasPhi.end(); ++it2 ) {
00373         if (((*it2)->globalPosition().phi() < seed_maxX[NNN] + dXclusBoxMax) && ((*it2)->globalPosition().phi() > seed_minX[NNN] - dXclusBoxMax)) {
00374           if(
00375              !( 
00376                (
00377                 ((*it2)->globalPosition().theta() + 0.3) < (seed_minY[NNN] - dYclusBoxMax) 
00378                 && 
00379                 ((*it2)->globalPosition().theta() - 0.3) < (seed_minY[NNN] - dYclusBoxMax)
00380                 )
00381                ||
00382                (
00383                 ((*it2)->globalPosition().theta() + 0.3) > (seed_maxY[NNN] + dYclusBoxMax) 
00384                 && 
00385                 ((*it2)->globalPosition().theta() - 0.3) > (seed_maxY[NNN] + dYclusBoxMax)
00386                 )
00387                )
00388              ) { // we have checked that the 2Dsegment is within tight phi boundaries and loose theta boundaries of the current cluster -> add it
00389             seeds[NNN].push_back((*it2)); // warning - neeed eta/theta switch here
00390             
00391           }
00392         }
00393       }
00394     } // DT2D_hastPhi loop
00395     
00396     // put RPC loop here
00397     int secondCh = 0;
00398     DetId detId_prev;
00399     if(seeds[NNN].size()>1){// actually we should check how many chambers with measurements are present
00400       for(unsigned int iRH = 0 ;iRH<seeds[NNN].size() ;++iRH){
00401         if( iRH && detId_prev != seeds[NNN][iRH]->hit()->geographicalId()){
00402           ++secondCh;
00403           break;
00404         }
00405         detId_prev = seeds[NNN][iRH]->hit()->geographicalId();
00406       }
00407     }
00408 
00409     if (useRPCs && !secondCh && !tooCloseClusters) {
00410       for (MuonRecHitContainer::const_iterator it2 = muonRecHits_RPC.begin(); it2 != muonRecHits_RPC.end(); ++it2 ) {
00411         if (((*it2)->globalPosition().phi() < seed_maxX[NNN] + dXclusBoxMax) && ((*it2)->globalPosition().phi() > seed_minX[NNN] - dXclusBoxMax)) {
00412           if(
00413              !( 
00414                (
00415                 ((*it2)->globalPosition().theta() + 0.3) < (seed_minY[NNN] - dYclusBoxMax) 
00416                 && 
00417                 ((*it2)->globalPosition().theta() - 0.3) < (seed_minY[NNN] - dYclusBoxMax)
00418                 )
00419                ||
00420                (
00421                 ((*it2)->globalPosition().theta() + 0.3) > (seed_maxY[NNN] + dYclusBoxMax) 
00422                 && 
00423                 ((*it2)->globalPosition().theta() - 0.3) > (seed_maxY[NNN] + dYclusBoxMax)
00424                 )
00425                )
00426              ) { // we have checked that the 2Dsegment is within tight phi boundaries and loose theta boundaries of the current cluster -> add it
00427             seeds[NNN].push_back((*it2)); // warning - neeed eta/theta switch here
00428             
00429           }
00430         }
00431       }
00432     } // RPC loop
00433   }
00434   
00435   // hand over the final seeds to the output
00436   // would be more elegant if we could do the above step with 
00437   // erasing the merged ones, rather than the 
00438   for(unsigned int NNN = 0; NNN < seeds.size(); ++NNN) {
00439     if(running_meanX[NNN] == 999999.) continue; //skip seeds that have been marked as used up in merging
00440     //std::cout<<"Next Cluster..."<<std::endl;
00441     segments_clusters.push_back(seeds[NNN]);
00442   }
00443 }
00444 
00445 
00446 bool SETPatternRecognition::segmentCleaning(const DetId & detId, 
00447                                             const LocalPoint& localPosition, const LocalError& localError,
00448                                             const LocalVector& localDirection, const LocalError& localDirectionError,
00449                                             const double& chi2, const int& ndf){
00450   // drop segments which are "bad"
00451   bool dropTheSegment = true;
00452   const GeomDet* geomDet = theService->trackingGeometry()->idToDet( detId );
00453   // only segments whithin the boundaries of the chamber
00454   bool insideCh = geomDet->surface().bounds().inside(localPosition, localError,outsideChamberErrorScale);
00455 
00456   // Don't use segments (nearly) parallel to the chamberi;
00457   // the direction vector is normalized (R=1)  
00458   bool parallelSegment = fabs(localDirection.z())>minLocalSegmentAngle? false: true;
00459 
00460   if(insideCh && !parallelSegment){
00461     dropTheSegment = false;
00462   }
00463   // use chi2 too? (DT, CSCs, RPCs; 2D, 4D;...)
00464 
00465 
00466   return dropTheSegment;
00467 }