CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/RecoLocalMuon/DTSegment/src/DTCombinatorialPatternReco.cc

Go to the documentation of this file.
00001 
00009 /* This Class Header */
00010 #include "RecoLocalMuon/DTSegment/src/DTCombinatorialPatternReco.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 
00027 
00028 /* C++ Headers */
00029 #include <iterator>
00030 using namespace std;
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032 
00033 /* ====================================================================== */
00034 
00036 DTCombinatorialPatternReco::DTCombinatorialPatternReco(const edm::ParameterSet& pset) : 
00037 DTRecSegment2DBaseAlgo(pset), theAlgoName("DTCombinatorialPatternReco")
00038 {
00039   theMaxAllowedHits = pset.getParameter<unsigned int>("MaxAllowedHits"); // 100
00040   theAlphaMaxTheta = pset.getParameter<double>("AlphaMaxTheta");// 0.1 ;
00041   theAlphaMaxPhi = pset.getParameter<double>("AlphaMaxPhi");// 1.0 ;
00042   debug = pset.getUntrackedParameter<bool>("debug"); //true;
00043   theUpdator = new DTSegmentUpdator(pset);
00044   theCleaner = new DTSegmentCleaner(pset);
00045   string theHitAlgoName = pset.getParameter<string>("recAlgo");
00046   usePairs = !(theHitAlgoName=="DTNoDriftAlgo");
00047 }
00048 
00050 DTCombinatorialPatternReco::~DTCombinatorialPatternReco() {
00051   delete theUpdator;
00052   delete theCleaner;
00053 }
00054 
00055 /* Operations */ 
00056 edm::OwnVector<DTSLRecSegment2D>
00057 DTCombinatorialPatternReco::reconstruct(const DTSuperLayer* sl,
00058                                         const std::vector<DTRecHit1DPair>& pairs){
00059 
00060   edm::OwnVector<DTSLRecSegment2D> result;
00061   vector<DTHitPairForFit*> hitsForFit = initHits(sl, pairs);
00062 
00063   vector<DTSegmentCand*> candidates = buildSegments(sl, hitsForFit);
00064 
00065   vector<DTSegmentCand*>::const_iterator cand=candidates.begin();
00066   while (cand<candidates.end()) {
00067     
00068     DTSLRecSegment2D *segment = (**cand);
00069 
00070     theUpdator->update(segment);
00071 
00072     result.push_back(segment);
00073 
00074     if (debug) {
00075       cout<<"Reconstructed 2D segments "<< result.back() <<endl;
00076     }
00077 
00078     delete *(cand++); // delete the candidate!
00079   }
00080 
00081   for (vector<DTHitPairForFit*>::iterator it = hitsForFit.begin(), ed = hitsForFit.end(); 
00082         it != ed; ++it) delete *it;
00083 
00084   return result;
00085 }
00086 
00087 void DTCombinatorialPatternReco::setES(const edm::EventSetup& setup){
00088   // Get the DT Geometry
00089   setup.get<MuonGeometryRecord>().get(theDTGeometry);
00090   theUpdator->setES(setup);
00091 }
00092 
00093 vector<DTHitPairForFit*>
00094 DTCombinatorialPatternReco::initHits(const DTSuperLayer* sl,
00095                                      const std::vector<DTRecHit1DPair>& hits){  
00096   
00097   vector<DTHitPairForFit*> result;
00098   for (vector<DTRecHit1DPair>::const_iterator hit=hits.begin();
00099        hit!=hits.end(); ++hit) {
00100     result.push_back(new DTHitPairForFit(*hit, *sl, theDTGeometry));
00101   }
00102   return result;
00103 }
00104 
00105 vector<DTSegmentCand*>
00106 DTCombinatorialPatternReco::buildSegments(const DTSuperLayer* sl,
00107                                           const std::vector<DTHitPairForFit*>& hits){
00108   // clear the patterns tried
00109   if (debug) {
00110       cout << "theTriedPattern.size is " << theTriedPattern.size() << "\n";
00111   }
00112   theTriedPattern.clear();
00113 
00114   typedef vector<DTHitPairForFit*> hitCont;
00115   typedef hitCont::const_iterator  hitIter;
00116   vector<DTSegmentCand*> result;
00117   
00118   if(debug) {
00119     cout << "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       //if ( (*lastHit)->id().layerId() == (*firstHit)->id().layerId() ) continue; // hits must be in different layers!
00143       // hits must nor in the same nor in adiacent layers
00144       if ( fabs((*lastHit)->id().layerId()-(*firstHit)->id().layerId())<=1 ) continue;
00145       if(debug) {
00146         cout << "Selected these two hits pair " << endl;
00147         cout << "First " << *(*firstHit) << " Layer Id: " << (*firstHit)->id().layerId() << endl;
00148         cout << "Last "  << *(*lastHit)  << " Layer Id: " << (*lastHit)->id().layerId()  << endl;
00149       }
00150 
00151       GlobalPoint IP;
00152       float DAlphaMax;
00153       if ((sl->id()).superlayer()==2)  // Theta SL
00154         DAlphaMax=theAlphaMaxTheta;
00155       else // Phi SL
00156         DAlphaMax=theAlphaMaxPhi;
00157 
00158       DTEnums::DTCellSide codes[2]={DTEnums::Right, DTEnums::Left};
00159       for (int firstLR=0; firstLR<2; ++firstLR) {
00160         for (int lastLR=0; lastLR<2; ++lastLR) {
00161           // TODO move the global transformation in the DTHitPairForFit class
00162           // when it will be moved I will able to remove the sl from the input parameter
00163           GlobalPoint gposFirst=sl->toGlobal( (*firstHit)->localPosition(codes[firstLR]) );
00164           GlobalPoint gposLast= sl->toGlobal( (*lastHit)->localPosition(codes[lastLR]) );
00165 
00166           GlobalVector gvec=gposLast-gposFirst;
00167           GlobalVector gvecIP=gposLast-IP;
00168 
00169           // difference in angle measured
00170           float DAlpha=fabs(gvec.theta()-gvecIP.theta());
00171 
00172           // cout << "DAlpha " << DAlpha << endl;
00173           if (DAlpha<DAlphaMax) {
00174 
00175             // create a segment hypotesis
00176             // I don't need a true segment, just direction and position
00177             LocalPoint posIni = (*firstHit)->localPosition(codes[firstLR]);
00178             LocalVector dirIni = 
00179               ((*lastHit)->localPosition(codes[lastLR])-posIni).unit();
00180 
00181             // search for other compatible hits, with or without the L/R solved
00182             vector<AssPoint> assHits = findCompatibleHits(posIni, dirIni, hits);
00183             if(debug) 
00184               cout << "compatible hits " << assHits.size() << endl;
00185 
00186             // get the best segment with these hits: it's just one! 
00187             // (is it correct?)
00188             DTSegmentCand* seg = buildBestSegment(assHits, sl);
00189 
00190             if (seg) {
00191               if(debug) 
00192                 cout << "segment " << *seg<< endl;
00193 
00194               // check if the chi2 and #hits are ok
00195               if (!seg->good()) {
00196                 delete seg;
00197               } else { 
00198 
00199                 // remove duplicated segments (I know, would be better to do it before the
00200                 // fit...)
00201                 if (checkDoubleCandidates(result,seg)) {
00202                   // add to the vector of hypotesis
00203                   result.push_back(seg);
00204                   if(debug) 
00205                     cout << "result is now " << result.size() << endl;
00206                 } else { // delete it!
00207                   delete seg;
00208                   if(debug) 
00209                     cout << "already existing" << endl;
00210                 }
00211               }
00212             }
00213           }
00214         }
00215       }
00216     }
00217   }
00218   if (debug) {
00219     for (vector<DTSegmentCand*>::const_iterator seg=result.begin();
00220          seg!=result.end(); ++seg) 
00221       cout << *(*seg) << endl;
00222   }
00223 
00224   // now I have a couple of segment hypotesis, should check for ghost
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   return result;
00234 }
00235 
00236 
00237 vector<DTCombinatorialPatternReco::AssPoint>
00238 DTCombinatorialPatternReco::findCompatibleHits(const LocalPoint& posIni,
00239                                                const LocalVector& dirIni,
00240                                                const vector<DTHitPairForFit*>& hits) {
00241   if (debug) cout << "Pos: " << posIni << " Dir: "<< dirIni << endl;
00242   vector<AssPoint> result;
00243 
00244   // counter to early-avoid double counting in hits pattern
00245   TriedPattern tried;
00246   int nCompatibleHits=0;
00247 
00248   typedef vector<DTHitPairForFit*> hitCont;
00249   typedef hitCont::const_iterator  hitIter;
00250   for (hitIter hit=hits.begin(); hit!=hits.end(); ++hit) {
00251     pair<bool,bool> isCompatible = (*hit)->isCompatible(posIni, dirIni);
00252     if (debug) 
00253       cout << "isCompatible " << isCompatible.first << " " <<
00254         isCompatible.second << endl;
00255 
00256     // if only one of the two is compatible, then the LR is assigned,
00257     // otherwise is undefined
00258 
00259     DTEnums::DTCellSide lrcode;
00260     if (isCompatible.first && isCompatible.second) {
00261       usePairs ? lrcode=DTEnums::undefLR : lrcode=DTEnums::Left ; // if not usePairs then only use single side 
00262       tried.push_back(3);
00263       nCompatibleHits++;
00264     }
00265     else if (isCompatible.first) {
00266       lrcode=DTEnums::Left;
00267       tried.push_back(2);
00268       nCompatibleHits++;
00269     }
00270     else if (isCompatible.second) {
00271       lrcode=DTEnums::Right;
00272       tried.push_back(1);
00273       nCompatibleHits++;
00274     }
00275     else {
00276       tried.push_back(0);
00277       continue; // neither is compatible
00278     }
00279     result.push_back(AssPoint(*hit, lrcode));
00280   }
00281 
00282 
00283   // check if too few associated hits or pattern already tried
00284   // TODO: two different if for nCompatibleHits < 3 =>printout and find ->
00285   // printour
00286   if (nCompatibleHits < 3) {
00287       if (debug) {
00288           cout << "Less than 3 hits " ;
00289           copy(tried.begin(), tried.end(), ostream_iterator<int>(std::cout));
00290           cout << endl;
00291       }
00292       // No need to insert into the list of patterns, as you will never search for it.
00293       //theTriedPattern.insert(tried);
00294   } else {
00295       // try to insert, return bool if already there
00296       bool isnew = theTriedPattern.insert(tried).second;
00297       if (isnew) {
00298           if (debug) {
00299               cout << "NOT Already tried " ;
00300               copy(tried.begin(), tried.end(), ostream_iterator<int>(std::cout));
00301               cout << endl;
00302           }
00303       } else {
00304           if (debug) {
00305               cout << "Already tried " ;
00306               copy(tried.begin(), tried.end(), ostream_iterator<int>(std::cout));
00307               cout << endl;
00308 
00309           }
00310           // empty the result vector
00311           result.clear();
00312       }
00313   }
00314   return result;
00315 }
00316 
00317 DTSegmentCand*
00318 DTCombinatorialPatternReco::buildBestSegment(std::vector<AssPoint>& hits,
00319                                              const DTSuperLayer* sl) {
00320   if (hits.size()<3) {
00321     //cout << "buildBestSegment: hits " << hits.size()<< endl;
00322     return 0; // a least 3 point
00323   }
00324 
00325   // hits with defined LR
00326   vector<AssPoint> points;
00327 
00328   // without: I store both L and R, a deque since I need front insertion and
00329   // deletion
00330   deque<DTHitPairForFit* > pointsNoLR; 
00331 
00332   // first add only the hits with LR assigned
00333   for (vector<AssPoint>::const_iterator hit=hits.begin();
00334        hit!=hits.end(); ++hit) {
00335     if ((*hit).second != DTEnums::undefLR) {
00336       points.push_back(*hit);
00337     } else { // then also for the undef'd one
00338       pointsNoLR.push_back((*hit).first);
00339     }
00340   }
00341 
00342   if(debug) {
00343     cout << "points " << points.size() << endl;
00344     cout << "pointsNoLR " << pointsNoLR.size() << endl;
00345   }
00346 
00347   // build all possible candidates using L/R ambiguity
00348   vector<DTSegmentCand*> candidates ;
00349 
00350   buildPointsCollection(points, pointsNoLR, candidates, sl);
00351 
00352   if(debug)
00353     cout << "candidates " << candidates.size() << endl;
00354 
00355   // so now I have build a given number of segments, I should find the best one,
00356   // by #hits and chi2.
00357   vector<DTSegmentCand*>::const_iterator bestCandIter = candidates.end();
00358   double minChi2=999999.;
00359   unsigned int maxNumHits=0;
00360   for (vector<DTSegmentCand*>::const_iterator iter=candidates.begin();
00361        iter!=candidates.end(); ++iter) {
00362     if ((*iter)->nHits()==maxNumHits && (*iter)->chi2()<minChi2) {
00363       minChi2=(*iter)->chi2();
00364       bestCandIter=iter;
00365     } else if ((*iter)->nHits()>maxNumHits) {
00366       maxNumHits=(*iter)->nHits();
00367       minChi2=(*iter)->chi2();
00368       bestCandIter=iter;
00369     }
00370   }
00371 
00372   // delete all candidates but the best one!
00373   for (vector<DTSegmentCand*>::iterator iter=candidates.begin();
00374        iter!=candidates.end(); ++iter) if (iter!=bestCandIter) delete *iter;
00375 
00376        // return the best candate if any
00377        if (bestCandIter != candidates.end()) {
00378          return (*bestCandIter);
00379        }
00380        return 0;
00381 }
00382 
00383 void
00384 DTCombinatorialPatternReco::buildPointsCollection(vector<AssPoint>& points, 
00385                                                   deque<DTHitPairForFit*>& pointsNoLR, 
00386                                                   vector<DTSegmentCand*>& candidates,
00387                                                   const DTSuperLayer* sl) {
00388 
00389   if(debug) {
00390     cout << "buildPointsCollection " << endl;
00391     cout << "points: " << points.size() << " NOLR: " << pointsNoLR.size()<< endl;
00392   }
00393   if (pointsNoLR.size()>0) { // still unassociated points!
00394     DTHitPairForFit* unassHit = pointsNoLR.front();
00395     // try with the right
00396     if(debug)
00397       cout << "Right hit" << endl;
00398     points.push_back(AssPoint(unassHit, DTEnums::Right));
00399     pointsNoLR.pop_front();
00400     buildPointsCollection(points, pointsNoLR, candidates, sl);
00401     pointsNoLR.push_front((unassHit));
00402     points.pop_back();
00403 
00404     // try with the left
00405     if(debug)
00406       cout << "Left hit" << endl;
00407     points.push_back(AssPoint(unassHit, DTEnums::Left));
00408     pointsNoLR.pop_front();
00409     buildPointsCollection(points, pointsNoLR, candidates, sl);
00410     pointsNoLR.push_front((unassHit));
00411     points.pop_back();
00412   } else { // all associated
00413 
00414     if(debug) {
00415       cout << "The Hits were" << endl;
00416       copy(points.begin(), points.end(),
00417            ostream_iterator<DTSegmentCand::AssPoint>(std::cout));
00418       cout << "----" << endl;
00419       cout << "All associated " << endl;
00420     }
00421     DTSegmentCand::AssPointCont pointsSet;
00422 
00423     // for (vector<AssPoint>::const_iterator point=points.begin();
00424     //      point!=points.end(); ++point) 
00425     pointsSet.insert(points.begin(),points.end());
00426 
00427     if(debug) {
00428       cout << "The Hits are" << endl;
00429       copy(pointsSet.begin(), pointsSet.end(),
00430            ostream_iterator<DTSegmentCand::AssPoint>(std::cout));
00431       cout << "----" << endl;
00432     }
00433 
00434     DTSegmentCand* newCand = new DTSegmentCand(pointsSet,sl);
00435     if (theUpdator->fit(newCand)) candidates.push_back(newCand);
00436     else delete newCand; // bad seg, too few hits
00437   }
00438 }
00439 
00440 bool
00441 DTCombinatorialPatternReco::checkDoubleCandidates(vector<DTSegmentCand*>& cands,
00442                                                   DTSegmentCand* seg) {
00443   for (vector<DTSegmentCand*>::iterator cand=cands.begin();
00444        cand!=cands.end(); ++cand) 
00445     if (*(*cand)==*seg) return false;
00446   return true;
00447 }