CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/RecoMuon/MuonSeedGenerator/src/MuonSeedOrcaPatternRecognition.cc

Go to the documentation of this file.
00001 
00014 #include "RecoMuon/MuonSeedGenerator/src/MuonSeedOrcaPatternRecognition.h"
00015 
00016 #include "DataFormats/CSCRecHit/interface/CSCRecHit2DCollection.h"
00017 #include "DataFormats/CSCRecHit/interface/CSCRecHit2D.h"
00018 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
00019 #include "DataFormats/DTRecHit/interface/DTRecSegment4D.h"
00020 #include "DataFormats/Math/interface/deltaPhi.h"
00021 #include "DataFormats/Common/interface/Handle.h"
00022 
00023 #include "RecoMuon/TransientTrackingRecHit/interface/MuonTransientTrackingRecHit.h"
00024 #include "RecoMuon/TrackingTools/interface/MuonPatternRecoDumper.h"
00025 
00026 
00027 // Geometry
00028 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00029 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00030 
00031 #include "RecoMuon/MeasurementDet/interface/MuonDetLayerMeasurements.h"
00032 #include "RecoMuon/DetLayers/interface/MuonDetLayerGeometry.h"
00033 #include "RecoMuon/Records/interface/MuonRecoGeometryRecord.h"
00034 
00035 // Framework
00036 #include "FWCore/Framework/interface/EventSetup.h"
00037 #include "FWCore/Framework/interface/Event.h"
00038 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00039 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00040 #include "FWCore/Framework/interface/ESHandle.h"
00041 
00042 // C++
00043 #include <vector>
00044 
00045 using namespace std;
00046 
00047     const std::string metname = "Muon|RecoMuon|MuonSeedOrcaPatternRecognition";
00048 
00049 // Constructor
00050 MuonSeedOrcaPatternRecognition::MuonSeedOrcaPatternRecognition(const edm::ParameterSet& pset)
00051 : MuonSeedVPatternRecognition(pset),
00052   theCrackEtas(pset.getParameter<std::vector<double> >("crackEtas")),
00053   theCrackWindow(pset.getParameter<double>("crackWindow"))
00054 {
00055 }
00056 
00057 
00058 // reconstruct muon's seeds
00059 void MuonSeedOrcaPatternRecognition::produce(const edm::Event& event, const edm::EventSetup& eSetup,
00060                                              std::vector<MuonRecHitContainer> & result)
00061 {
00062   // divide the RecHits by DetLayer, in order to fill the
00063   // RecHitContainer like it was in ORCA
00064   
00065   // Muon Geometry - DT, CSC and RPC 
00066   edm::ESHandle<MuonDetLayerGeometry> muonLayers;
00067   eSetup.get<MuonRecoGeometryRecord>().get(muonLayers);
00068 
00069   // get the DT layers
00070   vector<DetLayer*> dtLayers = muonLayers->allDTLayers();
00071 
00072   // get the CSC layers
00073   vector<DetLayer*> cscForwardLayers = muonLayers->forwardCSCLayers();
00074   vector<DetLayer*> cscBackwardLayers = muonLayers->backwardCSCLayers();
00075     
00076   // Backward (z<0) EndCap disk
00077   const DetLayer* ME4Bwd = cscBackwardLayers[4];
00078   const DetLayer* ME3Bwd = cscBackwardLayers[3];
00079   const DetLayer* ME2Bwd = cscBackwardLayers[2];
00080   const DetLayer* ME12Bwd = cscBackwardLayers[1];
00081   const DetLayer* ME11Bwd = cscBackwardLayers[0];
00082   
00083   // Forward (z>0) EndCap disk
00084   const DetLayer* ME11Fwd = cscForwardLayers[0];
00085   const DetLayer* ME12Fwd = cscForwardLayers[1];
00086   const DetLayer* ME2Fwd = cscForwardLayers[2];
00087   const DetLayer* ME3Fwd = cscForwardLayers[3];
00088   const DetLayer* ME4Fwd = cscForwardLayers[4];
00089      
00090   // barrel
00091   const DetLayer* MB4DL = dtLayers[3];
00092   const DetLayer* MB3DL = dtLayers[2];
00093   const DetLayer* MB2DL = dtLayers[1];
00094   const DetLayer* MB1DL = dtLayers[0];
00095   
00096   // instantiate the accessor
00097   // Don not use RPC for seeding
00098   MuonDetLayerMeasurements muonMeasurements(theDTRecSegmentLabel.label(),theCSCRecSegmentLabel,edm::InputTag(),
00099                                             enableDTMeasurement,enableCSCMeasurement,false);
00100   double barreldThetaCut = 0.2;
00101   // still lose good muons to a tighter cut
00102   double endcapdThetaCut = 1.0;
00103   MuonRecHitContainer list9 = filterSegments(muonMeasurements.recHits(MB4DL,event), barreldThetaCut);
00104   MuonRecHitContainer list6 = filterSegments(muonMeasurements.recHits(MB3DL,event), barreldThetaCut);
00105   MuonRecHitContainer list7 = filterSegments(muonMeasurements.recHits(MB2DL,event), barreldThetaCut);
00106   MuonRecHitContainer list8 = filterSegments(muonMeasurements.recHits(MB1DL,event), barreldThetaCut);
00107 
00108   dumpLayer("MB4 ", list9);
00109   dumpLayer("MB3 ", list6);
00110   dumpLayer("MB2 ", list7);
00111   dumpLayer("MB1 ", list8);
00112 
00113   bool* MB1 = zero(list8.size());
00114   bool* MB2 = zero(list7.size());
00115   bool* MB3 = zero(list6.size());
00116 
00117   endcapPatterns(filterSegments(muonMeasurements.recHits(ME11Bwd,event), endcapdThetaCut),
00118                  filterSegments(muonMeasurements.recHits(ME12Bwd,event), endcapdThetaCut),
00119                  filterSegments(muonMeasurements.recHits(ME2Bwd,event), endcapdThetaCut),
00120                  filterSegments(muonMeasurements.recHits(ME3Bwd,event), endcapdThetaCut),
00121                  filterSegments(muonMeasurements.recHits(ME4Bwd,event), endcapdThetaCut),
00122                  list8, list7, list6,
00123                  MB1, MB2, MB3, result);
00124 
00125   endcapPatterns(filterSegments(muonMeasurements.recHits(ME11Fwd,event), endcapdThetaCut),
00126                  filterSegments(muonMeasurements.recHits(ME12Fwd,event), endcapdThetaCut),
00127                  filterSegments(muonMeasurements.recHits(ME2Fwd,event), endcapdThetaCut),
00128                  filterSegments(muonMeasurements.recHits(ME3Fwd,event), endcapdThetaCut),
00129                  filterSegments(muonMeasurements.recHits(ME4Fwd,event), endcapdThetaCut),
00130                  list8, list7, list6,
00131                  MB1, MB2, MB3, result);
00132 
00133 
00134   // ----------    Barrel only
00135   
00136   unsigned int counter = 0;
00137   if ( list9.size() < 100 ) {   // +v
00138     for (MuonRecHitContainer::iterator iter=list9.begin(); iter!=list9.end(); iter++ ){
00139       MuonRecHitContainer seedSegments;
00140       seedSegments.push_back(*iter);
00141       complete(seedSegments, list6, MB3);
00142       complete(seedSegments, list7, MB2);
00143       complete(seedSegments, list8, MB1);
00144       if(check(seedSegments)) result.push_back(seedSegments);
00145     }
00146   }
00147 
00148 
00149   if ( list6.size() < 100 ) {   // +v
00150     for ( counter = 0; counter<list6.size(); counter++ ){
00151       if ( !MB3[counter] ) { 
00152         MuonRecHitContainer seedSegments;
00153         seedSegments.push_back(list6[counter]);
00154         complete(seedSegments, list7, MB2);
00155         complete(seedSegments, list8, MB1);
00156         complete(seedSegments, list9);
00157         if(check(seedSegments)) result.push_back(seedSegments);
00158       }
00159     }
00160   }
00161 
00162 
00163   if ( list7.size() < 100 ) {   // +v
00164     for ( counter = 0; counter<list7.size(); counter++ ){
00165       if ( !MB2[counter] ) { 
00166         MuonRecHitContainer seedSegments;
00167         seedSegments.push_back(list7[counter]);
00168         complete(seedSegments, list8, MB1);
00169         complete(seedSegments, list9);
00170         complete(seedSegments, list6, MB3);
00171         if (seedSegments.size()>1 || 
00172            (seedSegments.size()==1 && seedSegments[0]->dimension()==4) )
00173         {
00174           result.push_back(seedSegments);
00175         }
00176       }
00177     }
00178   }
00179 
00180 
00181   if ( list8.size() < 100 ) {   // +v
00182     for ( counter = 0; counter<list8.size(); counter++ ){
00183       if ( !MB1[counter] ) { 
00184         MuonRecHitContainer seedSegments;
00185         seedSegments.push_back(list8[counter]);
00186         complete(seedSegments, list9);
00187         complete(seedSegments, list6, MB3);
00188         complete(seedSegments, list7, MB2);
00189         if (seedSegments.size()>1 ||
00190            (seedSegments.size()==1 && seedSegments[0]->dimension()==4) )
00191         {
00192           result.push_back(seedSegments);
00193         }
00194       }
00195     }
00196   }
00197 
00198   if ( MB3 ) delete [] MB3;
00199   if ( MB2 ) delete [] MB2;
00200   if ( MB1 ) delete [] MB1;
00201 
00202   if(result.empty()) 
00203   {
00204     // be a little stricter with single segment seeds
00205     barreldThetaCut = 0.2;
00206     endcapdThetaCut = 0.2;
00207 
00208     MuonRecHitContainer all = muonMeasurements.recHits(ME4Bwd,event);
00209     MuonRecHitContainer tmp = filterSegments(muonMeasurements.recHits(ME3Bwd,event), endcapdThetaCut);
00210     copy(tmp.begin(),tmp.end(),back_inserter(all));
00211 
00212     tmp = filterSegments(muonMeasurements.recHits(ME2Bwd,event), endcapdThetaCut);
00213     copy(tmp.begin(),tmp.end(),back_inserter(all));
00214 
00215     tmp = filterSegments(muonMeasurements.recHits(ME12Bwd,event), endcapdThetaCut);
00216     copy(tmp.begin(),tmp.end(),back_inserter(all));
00217 
00218     tmp = filterSegments(muonMeasurements.recHits(ME11Bwd,event), endcapdThetaCut);
00219     copy(tmp.begin(),tmp.end(),back_inserter(all));
00220 
00221     tmp = filterSegments(muonMeasurements.recHits(ME11Fwd,event), endcapdThetaCut);
00222     copy(tmp.begin(),tmp.end(),back_inserter(all));
00223 
00224     tmp = filterSegments(muonMeasurements.recHits(ME12Fwd,event), endcapdThetaCut);
00225     copy(tmp.begin(),tmp.end(),back_inserter(all));
00226 
00227     tmp = filterSegments(muonMeasurements.recHits(ME2Fwd,event), endcapdThetaCut);
00228     copy(tmp.begin(),tmp.end(),back_inserter(all));
00229 
00230     tmp = filterSegments(muonMeasurements.recHits(ME3Fwd,event), endcapdThetaCut);
00231     copy(tmp.begin(),tmp.end(),back_inserter(all));
00232 
00233     tmp = filterSegments(muonMeasurements.recHits(ME4Fwd,event), endcapdThetaCut);
00234     copy(tmp.begin(),tmp.end(),back_inserter(all));
00235 
00236     tmp = filterSegments(muonMeasurements.recHits(MB4DL,event), barreldThetaCut);
00237     copy(tmp.begin(),tmp.end(),back_inserter(all));
00238 
00239     tmp = filterSegments(muonMeasurements.recHits(MB3DL,event), barreldThetaCut);
00240     copy(tmp.begin(),tmp.end(),back_inserter(all));
00241 
00242     tmp = filterSegments(muonMeasurements.recHits(MB2DL,event), barreldThetaCut);
00243     copy(tmp.begin(),tmp.end(),back_inserter(all));
00244 
00245     tmp = filterSegments(muonMeasurements.recHits(MB1DL,event), barreldThetaCut);
00246     copy(tmp.begin(),tmp.end(),back_inserter(all));
00247 
00248     LogTrace(metname)<<"Number of segments: "<<all.size();
00249 
00250     for(MuonRecHitContainer::const_iterator segmentItr = all.begin();
00251         segmentItr != all.end(); ++segmentItr)
00252     {
00253       MuonRecHitContainer singleSegmentContainer;
00254       singleSegmentContainer.push_back(*segmentItr);
00255       result.push_back(singleSegmentContainer);
00256     }
00257   }
00258 
00259 }
00260 
00261 
00262 bool * MuonSeedOrcaPatternRecognition::zero(unsigned listSize)
00263 {
00264   bool * result = 0;
00265   if (listSize) {
00266     result = new bool[listSize]; 
00267     for ( size_t i=0; i<listSize; i++ ) result[i]=false;
00268   }
00269   return result;
00270 }
00271 
00272 
00273 void MuonSeedOrcaPatternRecognition::endcapPatterns(
00274   const MuonRecHitContainer & me11, const MuonRecHitContainer & me12,
00275   const MuonRecHitContainer & me2,  const MuonRecHitContainer & me3,
00276   const MuonRecHitContainer & me4,  const  MuonRecHitContainer & mb1,
00277   const MuonRecHitContainer & mb2,  const  MuonRecHitContainer & mb3,
00278   bool * MB1, bool * MB2, bool * MB3,
00279   std::vector<MuonRecHitContainer> & result)
00280 {
00281   dumpLayer("ME4 ", me4);
00282   dumpLayer("ME3 ", me3);
00283   dumpLayer("ME2 ", me2);
00284   dumpLayer("ME12 ", me12);
00285   dumpLayer("ME11 ", me11);
00286 
00287   std::vector<MuonRecHitContainer> patterns;
00288   MuonRecHitContainer crackSegments;
00289   rememberCrackSegments(me11, crackSegments);
00290   rememberCrackSegments(me12, crackSegments);
00291   rememberCrackSegments(me2,  crackSegments);
00292   rememberCrackSegments(me3,  crackSegments);
00293   rememberCrackSegments(me4,  crackSegments);
00294 
00295 
00296   MuonRecHitContainer list24 = me4;
00297   MuonRecHitContainer list23 = me3;
00298 
00299   MuonRecHitContainer list12 = me2;
00300 
00301   MuonRecHitContainer list22 = me12;
00302   MuonRecHitContainer list21 = me11;
00303 
00304   MuonRecHitContainer list11 = list21;
00305   MuonRecHitContainer list5 = list22;
00306   MuonRecHitContainer list13 = list23;
00307   MuonRecHitContainer list4 = list24;
00308 
00309   if ( list21.size() == 0 )  {
00310     list11 = list22; list5 = list21;
00311   }
00312 
00313   if ( list24.size() < list23.size() && list24.size() > 0 )  {
00314     list13 = list24; list4 = list23;
00315   }
00316 
00317   if ( list23.size() == 0 )  {
00318     list13 = list24; list4 = list23;
00319   }
00320 
00321   MuonRecHitContainer list1 = list11;
00322   MuonRecHitContainer list2 = list12;
00323   MuonRecHitContainer list3 = list13;
00324 
00325 
00326   if ( list12.size() == 0 )  {
00327     list3 = list12;
00328     if ( list11.size() <= list13.size() && list11.size() > 0 ) {
00329       list1 = list11; list2 = list13;}
00330     else { list1 = list13; list2 = list11;}
00331   }
00332 
00333   if ( list13.size() == 0 )  {
00334     if ( list11.size() <= list12.size() && list11.size() > 0 ) {
00335       list1 = list11; list2 = list12;}
00336     else { list1 = list12; list2 = list11;}
00337   }
00338 
00339   if ( list12.size() != 0 &&  list13.size() != 0 )  {
00340     if ( list11.size()<=list12.size() && list11.size()<=list13.size() && list11.size()>0 ) {   // ME 1
00341       if ( list12.size() > list13.size() ) {
00342         list2 = list13; list3 = list12;}
00343     }
00344     else if ( list12.size() <= list13.size() ) {                                   //  start with ME 2
00345       list1 = list12;
00346       if ( list11.size() <= list13.size() && list11.size() > 0 ) {
00347         list2 = list11; list3 = list13;}
00348       else { list2 = list13; list3 = list11;}
00349     }
00350     else {                                                                         //  start with ME 3
00351       list1 = list13;
00352       if ( list11.size() <= list12.size() && list11.size() > 0 ) {
00353         list2 = list11; list3 = list12;}
00354       else { list2 = list12; list3 = list11;}
00355     }
00356   }
00357 
00358 
00359   bool* ME2 = zero(list2.size());
00360   bool* ME3 = zero(list3.size());
00361   bool* ME4 = zero(list4.size());
00362   bool* ME5 = zero(list5.size());
00363 
00364 
00365   // creates list of compatible track segments
00366 
00367   for (MuonRecHitContainer::iterator iter = list1.begin(); iter!=list1.end(); iter++ ){
00368     if ( (*iter)->recHits().size() < 4 && list3.size() > 0 ) continue; // 3p.tr-seg. are not so good for starting
00369     MuonRecHitContainer seedSegments;
00370     seedSegments.push_back(*iter);
00371     complete(seedSegments, list2, ME2);
00372     complete(seedSegments, list3, ME3);
00373     complete(seedSegments, list4, ME4);
00374     complete(seedSegments, list5, ME5);
00375     complete(seedSegments, mb3, MB3);
00376     complete(seedSegments, mb2, MB2);
00377     complete(seedSegments, mb1, MB1);
00378     if(check(seedSegments)) patterns.push_back(seedSegments);
00379   }
00380 
00381 
00382   unsigned int counter;
00383 
00384   for ( counter = 0; counter<list2.size(); counter++ ){
00385 
00386     if ( !ME2[counter] ) {
00387       MuonRecHitContainer seedSegments;
00388       seedSegments.push_back(list2[counter]);
00389       complete(seedSegments, list3, ME3);
00390       complete(seedSegments, list4, ME4);
00391       complete(seedSegments, list5, ME5);
00392       complete(seedSegments, mb3, MB3);
00393       complete(seedSegments, mb2, MB2);
00394       complete(seedSegments, mb1, MB1);
00395       if(check(seedSegments)) patterns.push_back(seedSegments);
00396     }
00397   }
00398 
00399 
00400   if ( list3.size() < 20 ) {   // +v
00401     for ( counter = 0; counter<list3.size(); counter++ ){
00402       if ( !ME3[counter] ) {
00403         MuonRecHitContainer seedSegments;
00404         seedSegments.push_back(list3[counter]);
00405         complete(seedSegments, list4, ME4);
00406         complete(seedSegments, list5, ME5);
00407         complete(seedSegments, mb3, MB3);
00408         complete(seedSegments, mb2, MB2);
00409         complete(seedSegments, mb1, MB1);
00410         if(check(seedSegments)) patterns.push_back(seedSegments);
00411       }
00412     }
00413   }
00414 
00415   if ( list4.size() < 20 ) {   // +v
00416     for ( counter = 0; counter<list4.size(); counter++ ){
00417       if ( !ME4[counter] ) {
00418         MuonRecHitContainer seedSegments;
00419         seedSegments.push_back(list4[counter]);
00420         complete(seedSegments, list5, ME5);
00421         complete(seedSegments, mb3, MB3);
00422         complete(seedSegments, mb2, MB2);
00423         complete(seedSegments, mb1, MB1);
00424         if(check(seedSegments)) patterns.push_back(seedSegments);
00425       }
00426     }
00427   }
00428 
00429   if ( ME5 ) delete [] ME5;
00430   if ( ME4 ) delete [] ME4;
00431   if ( ME3 ) delete [] ME3;
00432   if ( ME2 ) delete [] ME2;
00433 
00434   if(!patterns.empty())
00435   {
00436     result.insert(result.end(), patterns.begin(), patterns.end());
00437   }
00438   else
00439   {
00440     if(!crackSegments.empty())
00441     {
00442        // make some single-segment seeds
00443        for(MuonRecHitContainer::const_iterator crackSegmentItr = crackSegments.begin();
00444            crackSegmentItr != crackSegments.end(); ++crackSegmentItr)
00445        {
00446           MuonRecHitContainer singleSegmentPattern;
00447           singleSegmentPattern.push_back(*crackSegmentItr);
00448           result.push_back(singleSegmentPattern);
00449        }
00450     }
00451   }
00452 }
00453 
00454 
00455 
00456 void MuonSeedOrcaPatternRecognition::complete(MuonRecHitContainer& seedSegments,
00457                                  const MuonRecHitContainer &recHits, bool* used) const {
00458 
00459   MuonRecHitContainer good_rhit;
00460   MuonPatternRecoDumper theDumper;
00461   //+v get all rhits compatible with the seed on dEta/dPhi Glob.
00462 
00463   ConstMuonRecHitPointer first = seedSegments[0]; // first rechit of seed
00464 
00465   GlobalPoint ptg2 = first->globalPosition(); // its global pos +v
00466 
00467   for (unsigned nr = 0; nr < recHits.size(); ++nr ){
00468     MuonRecHitPointer recHit(recHits[nr]);
00469     GlobalPoint ptg1(recHit->globalPosition());
00470     float deta = fabs (ptg1.eta()-ptg2.eta());
00471     // Geom::Phi should keep it in the range [-pi, pi]
00472     float dphi = fabs( deltaPhi(ptg1.phi(), ptg2.phi()) );
00473     // be a little more lenient in cracks
00474     bool crack = isCrack(recHit) || isCrack(first);
00475     //float detaWindow = 0.3;
00476     float detaWindow = crack ? 0.25 : 0.2;
00477     if ( deta > detaWindow || dphi > .25 ) {
00478       continue;
00479     }   // +vvp!!!
00480 
00481     good_rhit.push_back(recHit);
00482     if (used) markAsUsed(nr, recHits, used);
00483   }  // recHits iter
00484 
00485   // select the best rhit among the compatible ones (based on Dphi Glob & Dir)
00486   MuonRecHitPointer best=bestMatch(first, good_rhit);
00487   if(best && best->isValid() ) seedSegments.push_back(best);
00488 }
00489 
00490 
00491 
00492 MuonSeedOrcaPatternRecognition::MuonRecHitPointer
00493 MuonSeedOrcaPatternRecognition::bestMatch(const ConstMuonRecHitPointer & first,
00494                                           MuonRecHitContainer & good_rhit) const
00495 {
00496   MuonRecHitPointer best = 0;
00497   if(good_rhit.size() == 1) return good_rhit[0];
00498   double bestDiscrim = 10000.;
00499   for (MuonRecHitContainer::iterator iter=good_rhit.begin(); 
00500        iter!=good_rhit.end(); iter++)
00501   {
00502     double discrim = discriminator(first, *iter);
00503     if(discrim < bestDiscrim)
00504     {
00505       bestDiscrim = discrim;
00506       best = *iter;
00507     }
00508   }
00509   return best;
00510 }
00511 
00512 
00513 double MuonSeedOrcaPatternRecognition::discriminator(const ConstMuonRecHitPointer & first, MuonRecHitPointer & other) const
00514 {
00515   GlobalPoint gp1= first->globalPosition();
00516   GlobalPoint gp2= other->globalPosition();
00517   GlobalVector gd1 = first->globalDirection();
00518   GlobalVector gd2 = other->globalDirection();
00519   if(first->isDT() || other->isDT()) {
00520     return fabs(deltaPhi(gd1.phi(), gd2.phi()));
00521   }
00522 
00523   // penalize those 3-hit segments
00524   int nhits = other->recHits().size();
00525   int penalty = std::max(nhits-2, 1);
00526   float dphig = deltaPhi(gp1.phi(), gp2.phi());
00527   // ME1A has slanted wires, so matching theta position doesn't work well.
00528   if(isME1A(first) || isME1A(other)) {
00529     return fabs(dphig/penalty);
00530   }
00531 
00532   float dthetag = gp1.theta()-gp2.theta();
00533   float dphid2 = fabs(deltaPhi(gd2.phi(), gp2.phi()));
00534   if (dphid2 > M_PI*.5) dphid2 = M_PI - dphid2;  //+v
00535   float dthetad2 = gp2.theta()-gd2.theta();
00536   // for CSC, make a big chi-squared of relevant variables
00537   // FIXME for 100 GeV mnd above muons, this doesn't perform as well as
00538   // previous methods.  needs investigation.
00539   float chisq =  ((dphig/0.02)*(dphig/0.02)
00540                 + (dthetag/0.003)*(dthetag/0.003)
00541                 + (dphid2/0.06)*(dphid2/0.06)
00542                 + (dthetad2/0.08)*(dthetad2/0.08)
00543                 );
00544   return chisq / penalty;
00545 }
00546 
00547 
00548 bool MuonSeedOrcaPatternRecognition::check(const MuonRecHitContainer & segments)
00549 {
00550   return (segments.size() > 1);
00551 }
00552 
00553 
00554 void MuonSeedOrcaPatternRecognition::markAsUsed(int nr, const MuonRecHitContainer &recHits, bool* used) const
00555 {
00556   used[nr] = true;
00557   // if it's ME1A with two other segments in the container, mark the ghosts as used, too.
00558   if(recHits[nr]->isCSC())
00559   {
00560     CSCDetId detId(recHits[nr]->geographicalId().rawId());
00561     if(detId.ring() == 4)
00562     {
00563       std::vector<unsigned> chamberHitNs;
00564       for(unsigned i = 0; i < recHits.size(); ++i)
00565       {
00566         if(recHits[i]->geographicalId().rawId() == detId.rawId())
00567         {
00568           chamberHitNs.push_back(i);
00569         }
00570       }
00571       if(chamberHitNs.size() == 3)
00572       {
00573         for(unsigned i = 0; i < 3; ++i)
00574         {
00575           used[chamberHitNs[i]] = true;
00576         }
00577       }
00578     }
00579   }
00580 }
00581 
00582 
00583 bool MuonSeedOrcaPatternRecognition::isCrack(const ConstMuonRecHitPointer & segment) const
00584 {
00585   bool result = false;
00586   double absEta = fabs(segment->globalPosition().eta());
00587   for(std::vector<double>::const_iterator crackItr = theCrackEtas.begin();
00588       crackItr != theCrackEtas.end(); ++crackItr)
00589   {
00590     if(fabs(absEta-*crackItr) < theCrackWindow) {
00591       result = true;
00592     }
00593   }
00594   return result;
00595 }
00596 
00597 
00598 void MuonSeedOrcaPatternRecognition::rememberCrackSegments(const MuonRecHitContainer & segments,
00599                                                            MuonRecHitContainer & crackSegments) const
00600 {
00601   for(MuonRecHitContainer::const_iterator segmentItr = segments.begin(); 
00602       segmentItr != segments.end(); ++segmentItr)
00603   {
00604     if((**segmentItr).hit()->dimension() == 4 && isCrack(*segmentItr)) 
00605     {
00606        crackSegments.push_back(*segmentItr);
00607     }
00608   }
00609 }
00610 
00611 
00612 
00613 void MuonSeedOrcaPatternRecognition::dumpLayer(const char * name, const MuonRecHitContainer & segments) const
00614 {
00615   MuonPatternRecoDumper theDumper;
00616 
00617   LogTrace(metname) << name << std::endl;
00618   for(MuonRecHitContainer::const_iterator segmentItr = segments.begin();
00619       segmentItr != segments.end(); ++segmentItr)
00620   {
00621     LogTrace(metname) << theDumper.dumpMuonId((**segmentItr).geographicalId());
00622   }
00623 }
00624 
00625 
00626 MuonSeedOrcaPatternRecognition::MuonRecHitContainer 
00627 MuonSeedOrcaPatternRecognition::filterSegments(const MuonRecHitContainer & segments, double dThetaCut) const
00628 {
00629 MuonPatternRecoDumper theDumper;
00630   MuonRecHitContainer result;
00631   for(MuonRecHitContainer::const_iterator segmentItr = segments.begin();
00632       segmentItr != segments.end(); ++segmentItr)
00633   {
00634     double dtheta = (*segmentItr)->globalDirection().theta() -  (*segmentItr)->globalPosition().theta();
00635     if((*segmentItr)->isDT())
00636     {
00637       // only apply the cut to 4D segments
00638       if((*segmentItr)->dimension() == 2 || fabs(dtheta) < dThetaCut)
00639       {
00640         result.push_back(*segmentItr);
00641       }
00642       else 
00643       {
00644         LogTrace(metname) << "Cutting segment " << theDumper.dumpMuonId((**segmentItr).geographicalId()) << " because dtheta = " << dtheta;
00645       }
00646 
00647     }
00648     else if((*segmentItr)->isCSC()) 
00649     {
00650       if(fabs(dtheta) < dThetaCut)
00651       {
00652         result.push_back(*segmentItr);
00653       }
00654       else 
00655       {
00656          LogTrace(metname) << "Cutting segment " << theDumper.dumpMuonId((**segmentItr).geographicalId()) << " because dtheta = " << dtheta;
00657       }
00658     }
00659   }
00660   filterOverlappingChambers(result);
00661   return result;
00662 }
00663 
00664 
00665 void MuonSeedOrcaPatternRecognition::filterOverlappingChambers(MuonRecHitContainer & segments) const
00666 {
00667   if(segments.empty()) return;
00668   MuonPatternRecoDumper theDumper;
00669   // need to optimize cuts
00670   double dphiCut = 0.05;
00671   double detaCut = 0.05;
00672   std::vector<unsigned> toKill;
00673   std::vector<unsigned> me1aOverlaps;
00674   // loop over all segment pairs to see if there are two that match up in eta and phi
00675   // but from different chambers
00676   unsigned nseg = segments.size();
00677   for(unsigned i = 0; i < nseg-1; ++i)
00678   {
00679     GlobalPoint pg1 = segments[i]->globalPosition();
00680     for(unsigned j = i+1; j < nseg; ++j)
00681     {
00682       GlobalPoint pg2 = segments[j]->globalPosition();
00683       if(segments[i]->geographicalId().rawId() != segments[j]->geographicalId().rawId()
00684          && fabs(deltaPhi(pg1.phi(), pg2.phi())) < dphiCut
00685          && fabs(pg1.eta()-pg2.eta()) < detaCut)
00686       {
00687         LogTrace(metname) << "OVERLAP " << theDumper.dumpMuonId(segments[i]->geographicalId()) << " " <<
00688                                    theDumper.dumpMuonId(segments[j]->geographicalId());
00689         // see which one is best
00690         toKill.push_back( (countHits(segments[i]) < countHits(segments[j])) ? i : j);
00691         if(isME1A(segments[i]))
00692         {
00693           me1aOverlaps.push_back(i);
00694           me1aOverlaps.push_back(j);
00695         }
00696       }
00697     }
00698   }
00699 
00700   if(toKill.empty()) return;
00701 
00702   // try to kill ghosts assigned to overlaps
00703   for(unsigned i = 0; i < me1aOverlaps.size(); ++i)
00704   {
00705     DetId detId(segments[me1aOverlaps[i]]->geographicalId());
00706     vector<unsigned> inSameChamber;
00707     for(unsigned j = 0; j < nseg; ++j)
00708     {
00709       if(i != j && segments[j]->geographicalId() == detId)
00710       {
00711         inSameChamber.push_back(j);
00712       }
00713     }
00714     if(inSameChamber.size() == 2)
00715     {
00716       toKill.push_back(inSameChamber[0]);
00717       toKill.push_back(inSameChamber[1]);
00718     }
00719   }
00720   // now kill the killable
00721   MuonRecHitContainer result;
00722   for(unsigned i = 0; i < nseg; ++i)
00723   {
00724     if(std::find(toKill.begin(), toKill.end(), i) == toKill.end())
00725     {
00726       result.push_back(segments[i]);
00727     }
00728    
00729   }
00730   segments.swap(result);
00731 }
00732 
00733 
00734 bool MuonSeedOrcaPatternRecognition::isME1A(const ConstMuonRecHitPointer & segment) const
00735 {
00736   return segment->isCSC() && CSCDetId(segment->geographicalId()).ring() == 4;
00737 }
00738 
00739 
00740 int MuonSeedOrcaPatternRecognition::countHits(const MuonRecHitPointer & segment) const {
00741   int count = 0;
00742   vector<TrackingRecHit*> components = (*segment).recHits();
00743   for(vector<TrackingRecHit*>::const_iterator component = components.begin();
00744       component != components.end(); ++component)
00745   {
00746     int componentSize = (**component).recHits().size();
00747     count += (componentSize == 0) ? 1 : componentSize;
00748   }
00749   return count;
00750 }
00751 
00752