CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/RecoLocalMuon/DTSegment/src/DTMeantimerPatternReco.cc

Go to the documentation of this file.
00001 
00010 /* This Class Header */
00011 #include "RecoLocalMuon/DTSegment/src/DTMeantimerPatternReco.h"
00012 
00013 /* Collaborating Class Header */
00014 #include "FWCore/Framework/interface/ESHandle.h"
00015 #include "FWCore/Framework/interface/EventSetup.h"
00016 #include "DataFormats/Common/interface/OwnVector.h"
00017 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00018 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00019 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00020 #include "Geometry/DTGeometry/interface/DTLayer.h"
00021 #include "Geometry/DTGeometry/interface/DTSuperLayer.h"
00022 #include "DataFormats/DTRecHit/interface/DTSLRecSegment2D.h"
00023 #include "RecoLocalMuon/DTSegment/src/DTSegmentUpdator.h"
00024 #include "RecoLocalMuon/DTSegment/src/DTSegmentCleaner.h"
00025 #include "RecoLocalMuon/DTSegment/src/DTHitPairForFit.h"
00026 #include "RecoLocalMuon/DTSegment/src/DTSegmentCand.h"
00027 
00028 /* C++ Headers */
00029 #include <iterator>
00030 using namespace std;
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032 
00033 /* ====================================================================== */
00034 
00035 
00037 DTMeantimerPatternReco::DTMeantimerPatternReco(const edm::ParameterSet& pset) : 
00038 DTRecSegment2DBaseAlgo(pset), theAlgoName("DTMeantimerPatternReco")
00039 {
00040 
00041   theMaxAllowedHits = pset.getParameter<unsigned int>("MaxAllowedHits"); // 100
00042   theAlphaMaxTheta = pset.getParameter<double>("AlphaMaxTheta");// 0.1 ;
00043   theAlphaMaxPhi = pset.getParameter<double>("AlphaMaxPhi");// 1.0 ;
00044   theMaxT0 = pset.getParameter<double>("MaxT0");
00045   theMinT0 = pset.getParameter<double>("MinT0");
00046   theMaxChi2 = pset.getParameter<double>("MaxChi2");// 8.0 ;
00047   debug = pset.getUntrackedParameter<bool>("debug");
00048   theUpdator = new DTSegmentUpdator(pset);
00049   theCleaner = new DTSegmentCleaner(pset);
00050 }
00051 
00053 DTMeantimerPatternReco::~DTMeantimerPatternReco() {
00054   delete theUpdator;
00055   delete theCleaner;
00056 }
00057 
00058 /* Operations */ 
00059 edm::OwnVector<DTSLRecSegment2D>
00060 DTMeantimerPatternReco::reconstruct(const DTSuperLayer* sl,
00061                                         const std::vector<DTRecHit1DPair>& pairs){
00062 
00063   edm::OwnVector<DTSLRecSegment2D> result;
00064   vector<DTHitPairForFit*> hitsForFit = initHits(sl, pairs);
00065 
00066   vector<DTSegmentCand*> candidates = buildSegments(sl, hitsForFit);
00067 
00068   vector<DTSegmentCand*>::const_iterator cand=candidates.begin();
00069   while (cand<candidates.end()) {
00070     
00071     DTSLRecSegment2D *segment = (**cand);
00072 
00073     theUpdator->update(segment);
00074 
00075     if (debug) 
00076       cout<<"Reconstructed 2D segments "<<*segment<<endl;
00077     result.push_back(segment);
00078 
00079     delete *(cand++); // delete the candidate!
00080   }
00081 
00082   return result;
00083 }
00084 
00085 void DTMeantimerPatternReco::setES(const edm::EventSetup& setup){
00086   // Get the DT Geometry
00087   setup.get<MuonGeometryRecord>().get(theDTGeometry);
00088   theUpdator->setES(setup);
00089 }
00090 
00091 vector<DTHitPairForFit*>
00092 DTMeantimerPatternReco::initHits(const DTSuperLayer* sl,
00093                                      const std::vector<DTRecHit1DPair>& hits){  
00094   
00095   vector<DTHitPairForFit*> result;
00096   for (vector<DTRecHit1DPair>::const_iterator hit=hits.begin();
00097        hit!=hits.end(); ++hit) {
00098     result.push_back(new DTHitPairForFit(*hit, *sl, theDTGeometry));
00099   }
00100   return result;
00101 }
00102 
00103 vector<DTSegmentCand*>
00104 DTMeantimerPatternReco::buildSegments(const DTSuperLayer* sl,
00105                                           const std::vector<DTHitPairForFit*>& hits){
00106 
00107   typedef vector<DTHitPairForFit*> hitCont;
00108   typedef hitCont::const_iterator  hitIter;
00109   vector<DTSegmentCand*> result;
00110   DTEnums::DTCellSide codes[2]={DTEnums::Left, DTEnums::Right};
00111 
00112   if(debug) {
00113     cout << "buildSegments: " << sl->id() << " nHits " << hits.size() << endl;
00114     for (hitIter hit=hits.begin(); hit!=hits.end(); ++hit) cout << **hit<< " wire: " << (*hit)->id() << " DigiTime: " << (*hit)->digiTime() << endl;
00115   }
00116 
00117   // 10-Mar-2004 SL
00118   // put a protection against heavily populated chambers, for which the segment
00119   // building could lead to infinite memory usage...
00120   if (hits.size() > theMaxAllowedHits ) {
00121     if(debug) {
00122       cout << "Warning: this SuperLayer " << sl->id() << " has too many hits : "
00123         << hits.size() << " max allowed is " << theMaxAllowedHits << endl;
00124       cout << "Skipping segment reconstruction... " << endl;
00125     }
00126     return result;
00127   }
00128 
00129   GlobalPoint IP;
00130   float DAlphaMax;
00131   if ((sl->id()).superlayer()==2)  // Theta SL
00132     DAlphaMax=theAlphaMaxTheta;
00133   else // Phi SL
00134     DAlphaMax=theAlphaMaxPhi;
00135   
00136   vector<AssPoint> usedHits;
00137 
00138   // get two hits in different layers and see if there are other hits
00139   //  compatible with them
00140   for (hitCont::const_iterator firstHit=hits.begin(); firstHit!=hits.end();
00141        ++firstHit) {
00142     for (hitCont::const_reverse_iterator lastHit=hits.rbegin(); 
00143          (*lastHit)!=(*firstHit); ++lastHit) {
00144       
00145       // a geometrical sensibility cut for the two hits
00146       if (!geometryFilter((*firstHit)->id(),(*lastHit)->id())) continue;
00147          
00148       // create a set of hits for the fit (only the hits between the two selected ones)
00149       hitCont hitsForFit;
00150       for (hitCont::const_iterator tmpHit=firstHit+1; (*tmpHit)!=(*lastHit); tmpHit++) 
00151         if ((geometryFilter((*tmpHit)->id(),(*lastHit)->id())) 
00152             && (geometryFilter((*tmpHit)->id(),(*firstHit)->id()))) hitsForFit.push_back(*tmpHit);
00153 
00154       for (int firstLR=0; firstLR<2; ++firstLR) {
00155         for (int lastLR=0; lastLR<2; ++lastLR) {
00156 
00157           // TODO move the global transformation in the DTHitPairForFit class
00158           // when it will be moved I will able to remove the sl from the input parameter
00159           GlobalPoint gposFirst=sl->toGlobal( (*firstHit)->localPosition(codes[firstLR]) );
00160           GlobalPoint gposLast= sl->toGlobal( (*lastHit)->localPosition(codes[lastLR]) );
00161           GlobalVector gvec=gposLast-gposFirst;
00162           GlobalVector gvecIP=gposLast-IP;
00163 
00164           // difference in angle measured
00165           float DAlpha=fabs(gvec.theta()-gvecIP.theta());
00166           if (DAlpha>DAlphaMax) continue;
00167           
00168           if(debug) {
00169               cout << "Selected hit pair:" << endl;
00170               cout << "  First " << *(*firstHit) << " Layer Id: " << (*firstHit)->id().layerId() << " Side: " << firstLR << " DigiTime: " << (*firstHit)->digiTime() << endl;
00171               cout << "  Last "  << *(*lastHit)  << " Layer Id: " << (*lastHit)->id().layerId()  << " Side: " << lastLR << " DigiTime: " << (*lastHit)->digiTime() <<  endl;
00172           }
00173         
00174           vector<AssPoint> assHits;
00175           // create a candidate hit list
00176           assHits.push_back(AssPoint(*firstHit,codes[firstLR]));
00177           assHits.push_back(AssPoint(*lastHit,codes[lastLR]));
00178 
00179           // run hit adding/segment building 
00180           maxfound = 3;
00181           addHits(sl,assHits,hitsForFit,result, usedHits);
00182         }
00183       }
00184     }
00185   }
00186 
00187   // now I have a couple of segment hypotheses, should check for ghost
00188   if (debug) {
00189     cout << "Result (before cleaning): " << result.size() << endl;
00190     for (vector<DTSegmentCand*>::const_iterator seg=result.begin(); seg!=result.end(); ++seg) 
00191       cout << *(*seg) << endl;
00192   }        
00193   result = theCleaner->clean(result);
00194   if (debug) {
00195     cout << "Result (after cleaning): " << result.size() << endl;
00196     for (vector<DTSegmentCand*>::const_iterator seg=result.begin();
00197          seg!=result.end(); ++seg) 
00198       cout << *(*seg) << endl;
00199   }
00200 
00201   return result;
00202 }
00203 
00204 
00205 
00206 void 
00207 DTMeantimerPatternReco::addHits(const DTSuperLayer* sl, vector<AssPoint>& assHits, const vector<DTHitPairForFit*>& hits, 
00208                                 vector<DTSegmentCand*> &result, vector<AssPoint>& usedHits) {
00209 
00210   typedef vector<DTHitPairForFit*> hitCont;
00211   double chi2l,chi2r,t0_corrl,t0_corrr;
00212   bool foundSomething = false;
00213 
00214   if (debug) 
00215     cout << "DTMeantimerPatternReco::addHit " << endl << "   Picked " << assHits.size() << " hits, " << hits.size() << " left." << endl;
00216   
00217   if (assHits.size()+hits.size()<maxfound) return;
00218           
00219   // loop over the remaining hits
00220   for (hitCont::const_iterator hit=hits.begin(); hit!=hits.end(); ++hit) {
00221     if (debug)
00222       cout << "     Trying B: " << **hit<< " wire: " << (*hit)->id() << endl;
00223             
00224     assHits.push_back(AssPoint(*hit, DTEnums::Left));
00225     bool left_ok=fitWithT0(assHits, chi2l, t0_corrl,0);
00226     assHits.pop_back();
00227 
00228     assHits.push_back(AssPoint(*hit, DTEnums::Right));
00229     bool right_ok=fitWithT0(assHits, chi2r, t0_corrr,0);
00230     assHits.pop_back();
00231 
00232     if (debug) {
00233       int nHits=assHits.size()+1;
00234       cout << "         Left:  t0_corr = " << t0_corrl << "ns  chi2/nHits = " << chi2l << "/" << nHits << "  ok: " << left_ok << endl;
00235       cout << "        Right:  t0_corr = " << t0_corrr << "ns  chi2/nHits = " << chi2r << "/" << nHits << "  ok: " << right_ok << endl;
00236     }
00237 
00238     if (!left_ok && !right_ok) continue;
00239 
00240     foundSomething = true;    
00241 
00242     // prepare the hit set for the next search, start from the other side                        
00243     hitCont hitsForFit;
00244     for (hitCont::const_iterator tmpHit=hit+1; tmpHit!=hits.end(); tmpHit++) 
00245       if (geometryFilter((*tmpHit)->id(),(*hit)->id())) hitsForFit.push_back(*tmpHit); 
00246       
00247     reverse(hitsForFit.begin(),hitsForFit.end());
00248 
00249     // choose only one - left or right
00250     if (assHits.size()>3 && left_ok && right_ok) {
00251       if (chi2l<chi2r-0.1) right_ok=false; else
00252         if (chi2r<chi2l-0.1) left_ok=false;
00253     }
00254     if (left_ok) { 
00255       assHits.push_back(AssPoint(*hit, DTEnums::Left));
00256       addHits(sl,assHits,hitsForFit,result,usedHits);
00257       assHits.pop_back();
00258     }
00259     if (right_ok) { 
00260       assHits.push_back(AssPoint(*hit, DTEnums::Right));
00261       addHits(sl,assHits,hitsForFit,result,usedHits);
00262       assHits.pop_back();
00263     }
00264   }
00265 
00266   if (foundSomething) return;
00267 
00268   // If we already have a segment with more hits from this hit pair - don't save this one.  
00269   if (assHits.size()<maxfound) return;
00270 
00271   // Check if semgent Ok, calculate chi2
00272   if (!fitWithT0(assHits, chi2l, t0_corrl,debug)) return;
00273 
00274   // If no more iterations - store the current segment
00275 
00276   DTSegmentCand::AssPointCont pointsSet;
00277   pointsSet.insert(assHits.begin(),assHits.end());
00278   DTSegmentCand* seg = new DTSegmentCand(pointsSet,sl);
00279   theUpdator->fit(seg);
00280 
00281   if (seg) {
00282     for (vector<AssPoint>::const_iterator hiti = assHits.begin()+1; hiti != assHits.end()-1; hiti++)
00283       usedHits.push_back(*hiti);
00284 
00285     if (assHits.size()>maxfound) maxfound = assHits.size();
00286     seg->setChi2(chi2l); // we need to set the chi^2 so that the cleaner can pick the best segments
00287     if (debug) cout << "Segment built: " << *seg<< endl;
00288     if (checkDoubleCandidates(result,seg)) {
00289       result.push_back(seg);
00290       if (debug) cout << "   Result is now " << result.size() << endl;
00291     } else {
00292       if (debug) cout << "   Exists - skipping" << endl;
00293       delete seg;
00294     }
00295   }
00296 }
00297 
00298 
00299 bool
00300 DTMeantimerPatternReco::geometryFilter( const DTWireId first, const DTWireId second ) const
00301 {
00302 //  return true;
00303 
00304   const int layerLowerCut[4]={0,-1,-2,-2};
00305   const int layerUpperCut[4]={0, 2, 2, 3};
00306 //  const int layerLowerCut[4]={0,-2,-4,-5};
00307 //  const int layerUpperCut[4]={0, 3, 4, 6};
00308 
00309   // deal only with hits that are in the same SL
00310   if (first.layerId().superlayerId().superLayer()!=second.layerId().superlayerId().superLayer()) 
00311     return true;
00312     
00313   int deltaLayer=abs(first.layerId().layer()-second.layerId().layer());
00314 
00315   // drop hits in the same layer
00316   if (!deltaLayer) return false;
00317 
00318   // protection against unexpected layer numbering
00319   if (deltaLayer>3) { 
00320     cout << "*** WARNING! DT Layer numbers differ by more than 3! for hits: " << endl;
00321     cout << "             " << first << endl;
00322     cout << "             " << second << endl;
00323     return false;
00324   }
00325 
00326   // accept only hits in cells "not too far away"
00327   int deltaWire=first.wire()-second.wire();
00328   if (second.layerId().layer()%2==0) deltaWire=-deltaWire; // yet another trick to get it right...
00329   if ((deltaWire<=layerLowerCut[deltaLayer]) || (deltaWire>=layerUpperCut[deltaLayer])) return false;
00330 
00331   return true;
00332 }
00333 
00334 
00335 bool
00336 DTMeantimerPatternReco::fitWithT0(const vector<AssPoint> &assHits, double &chi2, double &t0_corr, const bool fitdebug) 
00337 {
00338   typedef vector < pair<double,double> > hitCoord;
00339   double a,b,coordError,x,y;
00340   double sx=0,sy=0,sxy=0,sxx=0,ssx=0,ssy=0,s=0,ss=0;
00341   int leftHits=0,rightHits=0;
00342 
00343   if (assHits.size()<3) return false;
00344 
00345   // I'm assuming the single hit error is the same for all hits...
00346   coordError=((*(assHits.begin())).first)->localPositionError().xx();
00347 
00348   for (vector<AssPoint>::const_iterator assHit=assHits.begin(); assHit!=assHits.end(); ++assHit) {
00349     if (coordError!=((*(assHits.begin())).first)->localPositionError().xx()) 
00350       cout << "   DTMeantimerPatternReco: Warning! Hit errors different within one segment!" << endl;
00351 
00352     x=((*assHit).first)->localPosition((*assHit).second).z();
00353     y=((*assHit).first)->localPosition((*assHit).second).x();
00354 
00355     sx+=x;
00356     sy+=y;
00357     sxy+=x*y;
00358     sxx+=x*x;
00359     s++;
00360     
00361     if ((*assHit).second==DTEnums::Left) {
00362       leftHits++;
00363       ssx+=x;
00364       ssy+=y;
00365       ss++;
00366     } else {
00367       rightHits++;
00368       ssx-=x;
00369       ssy-=y;
00370       ss--;
00371     }  
00372   }      
00373                     
00374   if (fitdebug && debug)
00375     cout  << "   DTMeantimerPatternReco::fitWithT0   Left hits: " << leftHits << "  Right hits: " << rightHits << endl;
00376 
00377   if (leftHits && rightHits) {
00378 
00379     double delta = ss*ss*sxx+s*sx*sx+s*ssx*ssx-s*s*sxx-2*ss*sx*ssx;
00380     t0_corr=0.;
00381 
00382     if (delta) {
00383       a=(ssy*s*ssx+sxy*ss*ss+sy*sx*s-sy*ss*ssx-ssy*sx*ss-sxy*s*s)/delta;
00384       b=(ssx*sy*ssx+sxx*ssy*ss+sx*sxy*s-sxx*sy*s-ssx*sxy*ss-sx*ssy*ssx)/delta;
00385       t0_corr=(ssx*s*sxy+sxx*ss*sy+sx*sx*ssy-sxx*s*ssy-sx*ss*sxy-ssx*sx*sy)/delta;
00386     } else return false;
00387   } else {
00388     double d = s*sxx - sx*sx;
00389     b = (sxx*sy- sx*sxy)/ d;
00390     a = (s*sxy - sx*sy) / d;
00391     t0_corr=0;
00392   }
00393 
00394 // Calculate the chi^2 of the hits AFTER t0 correction
00395 
00396   double chi,chi2_not0;  
00397   chi2=0;
00398   chi2_not0=0;
00399   
00400   for (vector<AssPoint>::const_iterator assHit=assHits.begin(); assHit!=assHits.end(); ++assHit) {
00401     x=((*assHit).first)->localPosition((*assHit).second).z();
00402     y=((*assHit).first)->localPosition((*assHit).second).x();
00403 
00404     chi=y-a*x-b;
00405     chi2_not0+=chi*chi/coordError;
00406 
00407     if ((*assHit).second==DTEnums::Left) chi-=t0_corr;
00408       else chi+=t0_corr;
00409       
00410     chi2+=chi*chi/coordError;
00411   }
00412 
00413   // For 3-hit segments ignore timing information
00414   if (assHits.size()<4) {
00415     chi2=chi2_not0;
00416 //    if (chi2<theMaxChi2) return true;
00417     if (chi2<200.) return true;
00418       else return false;
00419   }
00420 
00421   t0_corr/=-0.00543; // convert drift distance to time
00422 
00423   if (debug && fitdebug) 
00424   {
00425     cout << "      t0_corr = " << t0_corr << "ns  chi2/nHits = " << chi2 << "/" << assHits.size() << endl;
00426     
00427     if (((chi2/(assHits.size()-2)<theMaxChi2) && (t0_corr<theMaxT0) && (t0_corr>theMinT0)) || (assHits.size()==4)) {
00428       cout << "      a  = " << a  << "       b  = " << b  << endl;
00429       for (vector<AssPoint>::const_iterator assHit=assHits.begin(); assHit!=assHits.end(); ++assHit) {
00430         x=((*assHit).first)->localPosition((*assHit).second).z();
00431         y=((*assHit).first)->localPosition((*assHit).second).x();
00432         cout << "   z= " << x << "   x= " << y;
00433         if ((*assHit).second==DTEnums::Left) cout << "   x_corr= " << y+t0_corr*0.00543;
00434                                        else  cout << "   x_corr= " << y-t0_corr*0.00543;
00435         cout << "   seg= " << a*x+b << endl;
00436       }
00437     }
00438   }
00439 
00440   if ((chi2/(assHits.size()-2)<theMaxChi2) && (t0_corr<theMaxT0) && (t0_corr>theMinT0)) return true;
00441     else return false;
00442 }
00443 
00444 void 
00445 DTMeantimerPatternReco::rawFit(double &a, double &b, const vector< pair<double,double> > &hits) {
00446 
00447   double s=0,sx=0,sy=0,x,y;
00448   double sxx=0,sxy=0;
00449 
00450   a=b=0;
00451   if (hits.size()==0) return;
00452     
00453   if (hits.size()==1) {
00454     b=(*(hits.begin())).second;
00455   } else {
00456     for (unsigned int i = 0; i != hits.size(); i++) {
00457       x=hits[i].first;
00458       y=hits[i].second;
00459       sy += y;
00460       sxy+= x*y;
00461       s += 1.;
00462       sx += x;
00463       sxx += x*x;
00464     }
00465     // protect against a vertical line???
00466 
00467     double d = s*sxx - sx*sx;
00468     b = (sxx*sy- sx*sxy)/ d;
00469     a = (s*sxy - sx*sy) / d;
00470   }
00471 }
00472 
00473 bool
00474 DTMeantimerPatternReco::checkDoubleCandidates(vector<DTSegmentCand*>& cands,
00475                                                   DTSegmentCand* seg) {
00476   for (vector<DTSegmentCand*>::iterator cand=cands.begin();
00477        cand!=cands.end(); ++cand) {
00478     if (*(*cand)==*seg) return false;
00479     if (((*cand)->nHits()>=seg->nHits()) && ((*cand)->chi2ndof()<seg->chi2ndof()))
00480       if ((*cand)->nSharedHitPairs(*seg)>int(seg->nHits()-2)) return false;
00481   }    
00482   return true;
00483 }