CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_6/src/RecoLocalMuon/DTSegment/src/DTCombinatorialExtendedPatternReco.cc

Go to the documentation of this file.
00001 
00009 /* This Class Header */
00010 #include "RecoLocalMuon/DTSegment/src/DTCombinatorialExtendedPatternReco.h"
00011 
00012 /* Collaborating Class Header */
00013 #include "FWCore/Framework/interface/ESHandle.h"
00014 #include "FWCore/Framework/interface/EventSetup.h"
00015 #include "DataFormats/Common/interface/OwnVector.h"
00016 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00017 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00018 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00019 #include "Geometry/DTGeometry/interface/DTLayer.h"
00020 #include "Geometry/DTGeometry/interface/DTSuperLayer.h"
00021 #include "DataFormats/DTRecHit/interface/DTSLRecSegment2D.h"
00022 #include "RecoLocalMuon/DTSegment/src/DTSegmentUpdator.h"
00023 #include "RecoLocalMuon/DTSegment/src/DTSegmentCleaner.h"
00024 #include "RecoLocalMuon/DTSegment/src/DTHitPairForFit.h"
00025 #include "RecoLocalMuon/DTSegment/src/DTSegmentCand.h"
00026 #include "RecoLocalMuon/DTSegment/src/DTSegmentExtendedCand.h"
00027 
00028 
00029 /* C++ Headers */
00030 #include <iterator>
00031 using namespace std;
00032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00033 
00034 /* ====================================================================== */
00035 
00037 DTCombinatorialExtendedPatternReco::DTCombinatorialExtendedPatternReco(const edm::ParameterSet& pset) : 
00038 DTRecSegment2DBaseAlgo(pset), theAlgoName("DTCombinatorialExtendedPatternReco")
00039 {
00040   theMaxAllowedHits = pset.getParameter<unsigned int>("MaxAllowedHits"); // 100
00041   theAlphaMaxTheta = pset.getParameter<double>("AlphaMaxTheta");// 0.1 ;
00042   theAlphaMaxPhi = pset.getParameter<double>("AlphaMaxPhi");// 1.0 ;
00043   debug = pset.getUntrackedParameter<bool>("debug"); //true;
00044   theUpdator = new DTSegmentUpdator(pset);
00045   theCleaner = new DTSegmentCleaner(pset);
00046   string theHitAlgoName = pset.getParameter<string>("recAlgo");
00047   usePairs = !(theHitAlgoName=="DTNoDriftAlgo");
00048 }
00049 
00051 DTCombinatorialExtendedPatternReco::~DTCombinatorialExtendedPatternReco() {
00052 }
00053 
00054 /* Operations */ 
00055 edm::OwnVector<DTSLRecSegment2D>
00056 DTCombinatorialExtendedPatternReco::reconstruct(const DTSuperLayer* sl,
00057                                                 const std::vector<DTRecHit1DPair>& pairs){
00058 
00059   if(debug) cout << "DTCombinatorialExtendedPatternReco::reconstruct" << endl;
00060   theTriedPattern.clear();
00061   edm::OwnVector<DTSLRecSegment2D> result;
00062   vector<DTHitPairForFit*> hitsForFit = initHits(sl, pairs);
00063 
00064   vector<DTSegmentCand*> candidates = buildSegments(sl, hitsForFit);
00065 
00066   vector<DTSegmentCand*>::const_iterator cand=candidates.begin();
00067   while (cand<candidates.end()) {
00068     
00069     DTSLRecSegment2D *segment = (**cand);
00070 
00071     theUpdator->update(segment);
00072 
00073     result.push_back(segment);
00074 
00075     if (debug) {
00076       cout<<"Reconstructed 2D extended segments "<< result.back() <<endl;
00077     }
00078 
00079     delete *(cand++); // delete the candidate!
00080   }
00081 
00082   for (vector<DTHitPairForFit*>::iterator it = hitsForFit.begin(), ed = hitsForFit.end(); 
00083         it != ed; ++it) delete *it;
00084 
00085   return result;
00086 }
00087 
00088 void DTCombinatorialExtendedPatternReco::setES(const edm::EventSetup& setup){
00089   // Get the DT Geometry
00090   setup.get<MuonGeometryRecord>().get(theDTGeometry);
00091   theUpdator->setES(setup);
00092 }
00093 
00094 void DTCombinatorialExtendedPatternReco::setClusters(vector<DTSLRecCluster> clusters) {
00095   theClusters = clusters;
00096 }
00097 
00098 vector<DTHitPairForFit*>
00099 DTCombinatorialExtendedPatternReco::initHits(const DTSuperLayer* sl,
00100                                      const std::vector<DTRecHit1DPair>& hits){  
00101   
00102   vector<DTHitPairForFit*> result;
00103   for (vector<DTRecHit1DPair>::const_iterator hit=hits.begin();
00104        hit!=hits.end(); ++hit) {
00105     result.push_back(new DTHitPairForFit(*hit, *sl, theDTGeometry));
00106   }
00107   return result;
00108 }
00109 
00110 vector<DTSegmentCand*>
00111 DTCombinatorialExtendedPatternReco::buildSegments(const DTSuperLayer* sl,
00112                                           const std::vector<DTHitPairForFit*>& hits){
00113 
00114   typedef vector<DTHitPairForFit*> hitCont;
00115   typedef hitCont::const_iterator  hitIter;
00116   vector<DTSegmentCand*> result;
00117   
00118   if(debug) {
00119     cout << "DTCombinatorialExtendedPatternReco::buildSegments: " << sl->id() << " nHits " << hits.size() << endl;
00120     for (vector<DTHitPairForFit*>::const_iterator hit=hits.begin();
00121          hit!=hits.end(); ++hit) cout << **hit<< endl;
00122   }
00123 
00124   // 10-Mar-2004 SL
00125   // put a protection against heavily populated chambers, for which the segment
00126   // building could lead to infinite memory usage...
00127   if (hits.size() > theMaxAllowedHits ) {
00128     if(debug) {
00129       cout << "Warning: this SuperLayer " << sl->id() << " has too many hits : "
00130         << hits.size() << " max allowed is " << theMaxAllowedHits << endl;
00131       cout << "Skipping segment reconstruction... " << endl;
00132     }
00133     return result;
00134   }
00135 
00137   //  compatible with them
00138   for (hitCont::const_iterator firstHit=hits.begin(); firstHit!=hits.end();
00139        ++firstHit) {
00140     for (hitCont::const_reverse_iterator lastHit=hits.rbegin(); 
00141          (*lastHit)!=(*firstHit); ++lastHit) {
00142       // hits must nor in the same nor in adiacent layers
00143       if ( fabs((*lastHit)->id().layerId()-(*firstHit)->id().layerId())<=1 ) continue;
00144       if(debug) {
00145         cout << "Selected these two hits pair " << endl;
00146         cout << "First " << *(*firstHit) << " Layer Id: " << (*firstHit)->id().layerId() << endl;
00147         cout << "Last "  << *(*lastHit)  << " Layer Id: " << (*lastHit)->id().layerId()  << endl;
00148       }
00149 
00150       GlobalPoint IP;
00151       float DAlphaMax;
00152       if ((sl->id()).superlayer()==2)  // Theta SL
00153         DAlphaMax=theAlphaMaxTheta;
00154       else // Phi SL
00155         DAlphaMax=theAlphaMaxPhi;
00156 
00157       DTEnums::DTCellSide codes[2]={DTEnums::Right, DTEnums::Left};
00158       for (int firstLR=0; firstLR<2; ++firstLR) {
00159         for (int lastLR=0; lastLR<2; ++lastLR) {
00160           // TODO move the global transformation in the DTHitPairForFit class
00161           // when it will be moved I will able to remove the sl from the input parameter
00162           GlobalPoint gposFirst=sl->toGlobal( (*firstHit)->localPosition(codes[firstLR]) );
00163           GlobalPoint gposLast= sl->toGlobal( (*lastHit)->localPosition(codes[lastLR]) );
00164 
00165           GlobalVector gvec=gposLast-gposFirst;
00166           GlobalVector gvecIP=gposLast-IP;
00167 
00168           // difference in angle measured
00169           float DAlpha=fabs(gvec.theta()-gvecIP.theta());
00170 
00171           // cout << "DAlpha " << DAlpha << endl;
00172           if (DAlpha<DAlphaMax) {
00173 
00174             // create a segment hypotesis
00175             // I don't need a true segment, just direction and position
00176             LocalPoint posIni = (*firstHit)->localPosition(codes[firstLR]);
00177             LocalVector dirIni = 
00178               ((*lastHit)->localPosition(codes[lastLR])-posIni).unit();
00179 
00180             // search for other compatible hits, with or without the L/R solved
00181             vector<AssPoint> assHits = findCompatibleHits(posIni, dirIni, hits);
00182             if(debug) 
00183               cout << "compatible hits " << assHits.size() << endl;
00184 
00185             // here return an extended candidate (which _has_ the original
00186             // segment)
00187             DTSegmentExtendedCand* seg = buildBestSegment(assHits, sl);
00188 
00189             if (seg) {
00190               if(debug) 
00191                 cout << "segment " << *seg<< endl;
00192 
00193               // check if the chi2 and #hits are ok
00194               if (!seg->good()) { // good is reimplmented in extended segment
00195                 delete seg;
00196               } else { 
00197 
00198                 // remove duplicated segments 
00199                 if (checkDoubleCandidates(result,seg)) {
00200                   // add to the vector of hypotesis
00201                   // still work with extended segments
00202                   result.push_back(seg);
00203                   if(debug) 
00204                     cout << "result is now " << result.size() << endl;
00205                 } else { // delete it!
00206                   delete seg;
00207                   if(debug) 
00208                     cout << "already existing" << endl;
00209                 }
00210               }
00211             }
00212           }
00213         }
00214       }
00215     }
00216   }
00217   if (debug) {
00218     for (vector<DTSegmentCand*>::const_iterator seg=result.begin();
00219          seg!=result.end(); ++seg) 
00220       cout << *(*seg) << endl;
00221   }
00222 
00223   // now I have a couple of segment hypotesis, should check for ghost
00224   // still with extended candidates
00225   result = theCleaner->clean(result);
00226   if (debug) {
00227     cout << "result no ghost  " << result.size() << endl;
00228     for (vector<DTSegmentCand*>::const_iterator seg=result.begin();
00229          seg!=result.end(); ++seg) 
00230       cout << *(*seg) << endl;
00231   }
00232 
00233   // here, finally, I have to return the set of _original_ segments, not the
00234   // extended ones.
00235   return result;
00236 }
00237 
00238 
00239 vector<DTCombinatorialExtendedPatternReco::AssPoint>
00240 DTCombinatorialExtendedPatternReco::findCompatibleHits(const LocalPoint& posIni,
00241                                                const LocalVector& dirIni,
00242                                                const vector<DTHitPairForFit*>& hits) {
00243   if (debug) cout << "Pos: " << posIni << " Dir: "<< dirIni << endl;
00244   vector<AssPoint> result;
00245 
00246   // counter to early-avoid double counting in hits pattern
00247   vector<int> tried;
00248   int nCompatibleHits=0;
00249 
00250   typedef vector<DTHitPairForFit*> hitCont;
00251   typedef hitCont::const_iterator  hitIter;
00252   for (hitIter hit=hits.begin(); hit!=hits.end(); ++hit) {
00253     pair<bool,bool> isCompatible = (*hit)->isCompatible(posIni, dirIni);
00254     if (debug) 
00255       cout << "isCompatible " << isCompatible.first << " " <<
00256         isCompatible.second << endl;
00257 
00258     // if only one of the two is compatible, then the LR is assigned,
00259     // otherwise is undefined
00260 
00261     DTEnums::DTCellSide lrcode;
00262     if (isCompatible.first && isCompatible.second) {
00263       usePairs ? lrcode=DTEnums::undefLR : lrcode=DTEnums::Left ; // if not usePairs then only use single side 
00264       tried.push_back(3);
00265       nCompatibleHits++;
00266     }
00267     else if (isCompatible.first) {
00268       lrcode=DTEnums::Left;
00269       tried.push_back(2);
00270       nCompatibleHits++;
00271     }
00272     else if (isCompatible.second) {
00273       lrcode=DTEnums::Right;
00274       tried.push_back(1);
00275       nCompatibleHits++;
00276     }
00277     else {
00278       tried.push_back(0);
00279       continue; // neither is compatible
00280     }
00281     result.push_back(AssPoint(*hit, lrcode));
00282   }
00283   
00284 
00285   // check if too few associated hits or pattern already tried
00286   if ( nCompatibleHits < 3 || find(theTriedPattern.begin(), theTriedPattern.end(),tried) == theTriedPattern.end()) {
00287     theTriedPattern.push_back(tried);
00288   } else {
00289     if (debug) {
00290       vector<vector<int> >::const_iterator t=find(theTriedPattern.begin(),
00291                                                   theTriedPattern.end(),
00292                                                   tried);
00293       cout << "Already tried";
00294       copy((*t).begin(), (*t).end(), ostream_iterator<int>(std::cout));
00295       cout << endl;
00296     }
00297     // empty the result vector
00298     result.clear();
00299   }
00300   return result;
00301 }
00302 
00303 DTSegmentExtendedCand*
00304 DTCombinatorialExtendedPatternReco::buildBestSegment(std::vector<AssPoint>& hits,
00305                                                      const DTSuperLayer* sl) {
00306   if (debug) cout << "DTCombinatorialExtendedPatternReco::buildBestSegment " <<
00307     hits.size()  << endl;
00308   if (hits.size()<3) {
00309     //cout << "buildBestSegment: hits " << hits.size()<< endl;
00310     return 0; // a least 3 point
00311   }
00312 
00313   // hits with defined LR
00314   vector<AssPoint> points;
00315 
00316   // without: I store both L and R, a deque since I need front insertion and
00317   // deletion
00318   deque<DTHitPairForFit* > pointsNoLR; 
00319 
00320   // first add only the hits with LR assigned
00321   for (vector<AssPoint>::const_iterator hit=hits.begin();
00322        hit!=hits.end(); ++hit) {
00323     if ((*hit).second != DTEnums::undefLR) {
00324       points.push_back(*hit);
00325     } else { // then also for the undef'd one
00326       pointsNoLR.push_back((*hit).first);
00327     }
00328   }
00329 
00330   if(debug) {
00331     cout << "points " << points.size() << endl;
00332     cout << "pointsNoLR " << pointsNoLR.size() << endl;
00333   }
00334 
00335   // build all possible candidates using L/R ambiguity
00336   vector<DTSegmentCand*> candidates ;
00337 
00338   buildPointsCollection(points, pointsNoLR, candidates, sl);
00339 
00340   // here I try to add the external clusters and build a set of "extended
00341   // segment candidate
00342   vector<DTSegmentExtendedCand*> extendedCands = extendCandidates(candidates,
00343                                                                   sl); 
00344   if (debug) cout << "extended candidates " << extendedCands.size() << endl;
00345 
00346   // so now I have build a given number of segments, I should find the best one,
00347   // by #hits and chi2.
00348   vector<DTSegmentExtendedCand*>::const_iterator bestCandIter = extendedCands.end();
00349   double minChi2=999999.;
00350   unsigned int maxNumHits=0;
00351   for (vector<DTSegmentExtendedCand*>::const_iterator iter=extendedCands.begin();
00352        iter!=extendedCands.end(); ++iter) {
00353     if ((*iter)->nHits()==maxNumHits && (*iter)->chi2()<minChi2) {
00354       minChi2=(*iter)->chi2();
00355       bestCandIter=iter;
00356     } else if ((*iter)->nHits()>maxNumHits) {
00357       maxNumHits=(*iter)->nHits();
00358       minChi2=(*iter)->chi2();
00359       bestCandIter=iter;
00360     }
00361   }
00362 
00363   // delete all candidates but the best one!
00364   for (vector<DTSegmentExtendedCand*>::iterator iter=extendedCands.begin();
00365        iter!=extendedCands.end(); ++iter) 
00366     if (iter!=bestCandIter) delete (*iter);
00367 
00368   // return the best candate if any
00369   if (bestCandIter != extendedCands.end()) {
00370     return (*bestCandIter);
00371   }
00372   return 0;
00373 }
00374 
00375 void
00376 DTCombinatorialExtendedPatternReco::buildPointsCollection(vector<AssPoint>& points, 
00377                                                   deque<DTHitPairForFit*>& pointsNoLR, 
00378                                                   vector<DTSegmentCand*>& candidates,
00379                                                   const DTSuperLayer* sl) {
00380 
00381   if(debug) {
00382     cout << "DTCombinatorialExtendedPatternReco::buildPointsCollection " << endl;
00383     cout << "points: " << points.size() << " NOLR: " << pointsNoLR.size()<< endl;
00384   }
00385   if (pointsNoLR.size()>0) { // still unassociated points!
00386     DTHitPairForFit* unassHit = pointsNoLR.front();
00387     // try with the right
00388     if(debug)
00389       cout << "Right hit" << endl;
00390     points.push_back(AssPoint(unassHit, DTEnums::Right));
00391     pointsNoLR.pop_front();
00392     buildPointsCollection(points, pointsNoLR, candidates, sl);
00393     pointsNoLR.push_front((unassHit));
00394     points.pop_back();
00395 
00396     // try with the left
00397     if(debug)
00398       cout << "Left hit" << endl;
00399     points.push_back(AssPoint(unassHit, DTEnums::Left));
00400     pointsNoLR.pop_front();
00401     buildPointsCollection(points, pointsNoLR, candidates, sl);
00402     pointsNoLR.push_front((unassHit));
00403     points.pop_back();
00404   } else { // all associated
00405 
00406     if(debug) {
00407       cout << "The Hits were" << endl;
00408       copy(points.begin(), points.end(),
00409            ostream_iterator<DTSegmentCand::AssPoint>(std::cout));
00410       cout << "----" << endl;
00411       cout << "All associated " << endl;
00412     }
00413     DTSegmentCand::AssPointCont pointsSet;
00414 
00415     // for (vector<AssPoint>::const_iterator point=points.begin();
00416     //      point!=points.end(); ++point) 
00417     pointsSet.insert(points.begin(),points.end());
00418 
00419     if(debug) {
00420       cout << "The Hits are" << endl;
00421       copy(pointsSet.begin(), pointsSet.end(),
00422            ostream_iterator<DTSegmentCand::AssPoint>(std::cout));
00423       cout << "----" << endl;
00424     }
00425 
00426     DTSegmentCand* newCand = new DTSegmentCand(pointsSet,sl);
00427     if (theUpdator->fit(newCand)) candidates.push_back(newCand);
00428     else delete newCand; // bad seg, too few hits
00429   }
00430 }
00431 
00432 bool
00433 DTCombinatorialExtendedPatternReco::checkDoubleCandidates(vector<DTSegmentCand*>& cands,
00434                                                   DTSegmentCand* seg) {
00435   for (vector<DTSegmentCand*>::iterator cand=cands.begin();
00436        cand!=cands.end(); ++cand) 
00437     if (*(*cand)==*seg) return false;
00438   return true;
00439 }
00440 
00441 vector<DTSegmentExtendedCand*>
00442 DTCombinatorialExtendedPatternReco::extendCandidates(vector<DTSegmentCand*>& candidates,
00443                                                      const DTSuperLayer* sl) {
00444   if (debug) cout << "extendCandidates " << candidates.size() << endl;
00445   vector<DTSegmentExtendedCand*> result;
00446 
00447   // in case of phi SL just return
00448   if (sl->id().superLayer() != 2 ) {
00449     for (vector<DTSegmentCand*>:: const_iterator cand=candidates.begin();
00450        cand!=candidates.end(); ++cand) {
00451       DTSegmentExtendedCand* extendedCand = new DTSegmentExtendedCand(*cand);
00452       // and delete the original candidate
00453       delete *cand;
00454       result.push_back(extendedCand);
00455     }
00456     return result;
00457   }
00458 
00459   // first I have to select the cluster which are compatible with the actual
00460   // candidate, namely +/-1 sector/station/wheel 
00461   vector<DTSegmentExtendedCand::DTSLRecClusterForFit> clustersWithPos;
00462   if (debug) cout << "AllClustersWithPos " << theClusters.size() << endl;
00463   if(debug) cout << "SL:   " << sl->id() << endl;
00464   for (vector<DTSLRecCluster>::const_iterator clus=theClusters.begin();
00465        clus!=theClusters.end(); ++clus) {
00466     if(debug) cout << "CLUS: " << (*clus).superLayerId() << endl;
00467     if ((*clus).superLayerId().superLayer()==2 && closeSL(sl->id(),(*clus).superLayerId())) {
00468       // and then get their pos in the actual SL frame
00469       const DTSuperLayer* clusSl =
00470         theDTGeometry->superLayer((*clus).superLayerId());
00471       LocalPoint pos=sl->toLocal(clusSl->toGlobal((*clus).localPosition()));
00472       //LocalError err=sl->toLocal(clusSl->toGlobal((*clus).localPositionError()));
00473       LocalError err=(*clus).localPositionError();
00474       clustersWithPos.push_back(DTSegmentExtendedCand::DTSLRecClusterForFit(*clus, pos, err) );
00475     }
00476   }
00477   if (debug) cout << "closeClustersWithPos " << clustersWithPos.size() << endl;
00478 
00479   for (vector<DTSegmentCand*>:: const_iterator cand=candidates.begin();
00480        cand!=candidates.end(); ++cand) {
00481     // create an extended candidate
00482     DTSegmentExtendedCand* extendedCand = new DTSegmentExtendedCand(*cand);
00483     // and delete the original candidate
00484     delete *cand;
00485     // do this only for theta SL
00486     if (extendedCand->superLayer()->id().superLayer() == 2 ) {
00487       // first check compatibility between cand and clusForFit
00488       for (vector<DTSegmentExtendedCand::DTSLRecClusterForFit>::const_iterator
00489            exClus=clustersWithPos.begin(); exClus!=clustersWithPos.end(); ++exClus) {
00490         if (extendedCand->isCompatible(*exClus)) {
00491           if (debug) cout << "is compatible " << endl;
00492           // add compatible cluster
00493           extendedCand->addClus(*exClus);
00494         }
00495       }
00496       // fit the segment
00497       if (debug) cout << "extended cands nHits: " << extendedCand->nHits() <<endl;
00498       if (theUpdator->fit(extendedCand)) {
00499         // add to result
00500         result.push_back(extendedCand);
00501       } else {
00502         cout << "Bad fit" << endl;
00503         delete extendedCand;
00504       }
00505     } else { // Phi SuperLayer
00506       result.push_back(extendedCand);
00507     }
00508   }
00509 
00510   return result;
00511 }
00512 
00513 bool DTCombinatorialExtendedPatternReco::closeSL(const DTSuperLayerId& id1,
00514                                                  const DTSuperLayerId& id2) {
00515   if (id1==id2) return false;
00516   if (abs(id1.wheel()-id2.wheel())>1 ) return false;
00517   // take into account also sector 13 and 14
00518   int sec1 = ( id1.sector()==13 ) ? 4: id1.sector();
00519   sec1=(sec1==14)? 10: sec1;
00520   int sec2 = ( id2.sector()==13 ) ? 4: id2.sector();
00521   sec2=(sec2==14)? 10: sec2;
00522   // take into account also sector 1/12
00523   if (abs(sec1-sec2)>1 && abs(sec1-sec2)!=11 ) return false;
00524   //if (abs(id1.station()-id2.station())>1 ) return false;
00525   return true;
00526 }