CMS 3D CMS Logo

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

Generated on Tue Jun 9 17:43:53 2009 for CMSSW by  doxygen 1.5.4