CMS 3D CMS Logo

DTMeantimerPatternReco Class Reference

Algo for reconstructing 2d segment in DT using a combinatorial approach with a T0 estimation produced along the way. More...

#include <RecoLocalMuon/DTSegment/src/DTMeantimerPatternReco.h>

Inheritance diagram for DTMeantimerPatternReco:

DTRecSegment2DBaseAlgo

List of all members.

Public Member Functions

virtual std::string algoName () const
 return the algo name
 DTMeantimerPatternReco (const edm::ParameterSet &pset)
 Constructor.
virtual edm::OwnVector
< DTSLRecSegment2D
reconstruct (const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits)
 this function is called in the producer
virtual void setES (const edm::EventSetup &setup)
 Through this function the EventSetup is percolated to the objs which request it.
virtual ~DTMeantimerPatternReco ()
 Destructor.

Private Types

typedef std::pair
< DTHitPairForFit *,
DTEnums::DTCellSide
AssPoint

Private Member Functions

void addHits (const DTSuperLayer *sl, std::vector< AssPoint > &assHits, const std::vector< DTHitPairForFit * > &hits, std::vector< DTSegmentCand * > &result, std::vector< AssPoint > &usedHits)
std::vector< DTSegmentCand * > buildSegments (const DTSuperLayer *sl, const std::vector< DTHitPairForFit * > &hits)
bool checkDoubleCandidates (std::vector< DTSegmentCand * > &segs, DTSegmentCand *seg)
bool fitWithT0 (const std::vector< AssPoint > &assHits, double &chi2, double &t0_corr, const bool fitdebug)
bool geometryFilter (const DTWireId first, const DTWireId second) const
std::vector< DTHitPairForFit * > initHits (const DTSuperLayer *sl, const std::vector< DTRecHit1DPair > &hits)
void rawFit (double &a, double &b, const std::vector< std::pair< double, double > > &hits)

Private Attributes

bool debug
unsigned int maxfound
std::string theAlgoName
double theAlphaMaxPhi
double theAlphaMaxTheta
DTSegmentCleanertheCleaner
edm::ESHandle< DTGeometrytheDTGeometry
unsigned int theMaxAllowedHits
double theMaxChi2
double theMaxT0
double theMinT0
DTSegmentUpdatortheUpdator

Friends

class DTMeantimerPatternReco4D


Detailed Description

Algo for reconstructing 2d segment in DT using a combinatorial approach with a T0 estimation produced along the way.

Date
2008/03/10 11:18:20
Revision
1.2
Author:
Stefano Lacaprara - INFN Legnaro <stefano.lacaprara@pd.infn.it>

Riccardo Bellan - INFN TO <riccardo.bellan@cern.ch>

Piotr Traczyk - SINS Warsaw <ptraczyk@fuw.edu.pl>

Definition at line 43 of file DTMeantimerPatternReco.h.


Member Typedef Documentation

typedef std::pair<DTHitPairForFit*, DTEnums::DTCellSide> DTMeantimerPatternReco::AssPoint [private]

Definition at line 72 of file DTMeantimerPatternReco.h.


Constructor & Destructor Documentation

DTMeantimerPatternReco::DTMeantimerPatternReco ( const edm::ParameterSet pset  ) 

Constructor.

Definition at line 37 of file DTMeantimerPatternReco.cc.

References debug, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), theAlphaMaxPhi, theAlphaMaxTheta, theCleaner, theMaxAllowedHits, theMaxChi2, theMaxT0, theMinT0, and theUpdator.

00037                                                                           : 
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 }

DTMeantimerPatternReco::~DTMeantimerPatternReco (  )  [virtual]

Destructor.

Definition at line 53 of file DTMeantimerPatternReco.cc.

References theCleaner, and theUpdator.

00053                                                 {
00054   delete theUpdator;
00055   delete theCleaner;
00056 }


Member Function Documentation

void DTMeantimerPatternReco::addHits ( const DTSuperLayer sl,
std::vector< AssPoint > &  assHits,
const std::vector< DTHitPairForFit * > &  hits,
std::vector< DTSegmentCand * > &  result,
std::vector< AssPoint > &  usedHits 
) [private]

Referenced by buildSegments().

virtual std::string DTMeantimerPatternReco::algoName ( void   )  const [inline, virtual]

return the algo name

Implements DTRecSegment2DBaseAlgo.

Definition at line 61 of file DTMeantimerPatternReco.h.

References theAlgoName.

00061 { return theAlgoName; }

vector< DTSegmentCand * > DTMeantimerPatternReco::buildSegments ( const DTSuperLayer sl,
const std::vector< DTHitPairForFit * > &  hits 
) [private]

Definition at line 104 of file DTMeantimerPatternReco.cc.

References addHits(), DTSegmentCleaner::clean(), GenMuonPlsPt100GeV_cfg::cout, debug, lat::endl(), geometryFilter(), DTSuperLayer::id(), DTEnums::Left, maxfound, HLT_VtxMuL3::result, DTEnums::Right, seg, theAlphaMaxPhi, theAlphaMaxTheta, theCleaner, theMaxAllowedHits, PV3DBase< T, PVType, FrameType >::theta(), and GeomDet::toGlobal().

Referenced by DTMeantimerPatternReco4D::buildPhiSuperSegmentsCandidates(), and reconstruct().

00105                                                                                   {
00106 
00107   typedef vector<DTHitPairForFit*> hitCont;
00108   typedef hitCont::const_iterator  hitIter;
00109   vector<DTSegmentCand*> result;
00110   DTEnums::DTCellSide codes[2]={DTEnums::Right, DTEnums::Left};
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() << 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 << endl;
00171               cout << "  Last "  << *(*lastHit)  << " Layer Id: " << (*lastHit)->id().layerId()  << " Side: " << lastLR << 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 }

bool DTMeantimerPatternReco::checkDoubleCandidates ( std::vector< DTSegmentCand * > &  segs,
DTSegmentCand seg 
) [private]

bool DTMeantimerPatternReco::fitWithT0 ( const std::vector< AssPoint > &  assHits,
double &  chi2,
double &  t0_corr,
const bool  fitdebug 
) [private]

bool DTMeantimerPatternReco::geometryFilter ( const DTWireId  first,
const DTWireId  second 
) const [private]

Definition at line 301 of file DTMeantimerPatternReco.cc.

References funct::abs(), GenMuonPlsPt100GeV_cfg::cout, DTLayerId::layer(), DTWireId::layerId(), DTSuperLayerId::superLayer(), DTLayerId::superlayerId(), and DTWireId::wire().

Referenced by buildSegments().

00302 {
00303 //  return true;
00304 
00305   const int layerLowerCut[4]={0,-1,-2,-2};
00306   const int layerUpperCut[4]={0, 2, 2, 3};
00307 //  const int layerLowerCut[4]={0,-2,-4,-5};
00308 //  const int layerUpperCut[4]={0, 3, 4, 6};
00309 
00310   // deal only with hits that are in the same SL
00311   if (first.layerId().superlayerId().superLayer()!=second.layerId().superlayerId().superLayer()) 
00312     return true;
00313     
00314   int deltaLayer=abs(first.layerId().layer()-second.layerId().layer());
00315 
00316   // drop hits in the same layer
00317   if (!deltaLayer) return false;
00318 
00319   // protection against unexpected layer numbering
00320   if (deltaLayer>3) { 
00321     cout << "*** WARNING! DT Layer numbers differ by more than 3! for hits: " << endl;
00322     cout << "             " << first << endl;
00323     cout << "             " << second << endl;
00324     return false;
00325   }
00326 
00327   // accept only hits in cells "not too far away"
00328   int deltaWire=first.wire()-second.wire();
00329   if (second.layerId().layer()%2==0) deltaWire=-deltaWire; // yet another trick to get it right...
00330   if ((deltaWire<=layerLowerCut[deltaLayer]) || (deltaWire>=layerUpperCut[deltaLayer])) return false;
00331 
00332   return true;
00333 }

vector< DTHitPairForFit * > DTMeantimerPatternReco::initHits ( const DTSuperLayer sl,
const std::vector< DTRecHit1DPair > &  hits 
) [private]

Definition at line 92 of file DTMeantimerPatternReco.cc.

References HLT_VtxMuL3::result, and theDTGeometry.

Referenced by DTMeantimerPatternReco4D::buildPhiSuperSegmentsCandidates(), and reconstruct().

00093                                                                            {  
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 }

void DTMeantimerPatternReco::rawFit ( double &  a,
double &  b,
const std::vector< std::pair< double, double > > &  hits 
) [private]

edm::OwnVector< DTSLRecSegment2D > DTMeantimerPatternReco::reconstruct ( const DTSuperLayer sl,
const std::vector< DTRecHit1DPair > &  hits 
) [virtual]

this function is called in the producer

Implements DTRecSegment2DBaseAlgo.

Definition at line 60 of file DTMeantimerPatternReco.cc.

References buildSegments(), GenMuonPlsPt100GeV_cfg::cout, debug, lat::endl(), initHits(), edm::OwnVector< T, P >::push_back(), HLT_VtxMuL3::result, theUpdator, and DTSegmentUpdator::update().

Referenced by DTMeantimerPatternReco4D::reconstruct().

00061                                                                                {
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 }

void DTMeantimerPatternReco::setES ( const edm::EventSetup setup  )  [virtual]

Through this function the EventSetup is percolated to the objs which request it.

Implements DTRecSegment2DBaseAlgo.

Definition at line 85 of file DTMeantimerPatternReco.cc.

References edm::EventSetup::get(), DTSegmentUpdator::setES(), theDTGeometry, and theUpdator.

Referenced by DTMeantimerPatternReco4D::setES().

00085                                                             {
00086   // Get the DT Geometry
00087   setup.get<MuonGeometryRecord>().get(theDTGeometry);
00088   theUpdator->setES(setup);
00089 }


Friends And Related Function Documentation

friend class DTMeantimerPatternReco4D [friend]

Definition at line 70 of file DTMeantimerPatternReco.h.


Member Data Documentation

bool DTMeantimerPatternReco::debug [private]

Definition at line 110 of file DTMeantimerPatternReco.h.

Referenced by buildSegments(), DTMeantimerPatternReco(), and reconstruct().

unsigned int DTMeantimerPatternReco::maxfound [private]

Definition at line 114 of file DTMeantimerPatternReco.h.

Referenced by buildSegments().

std::string DTMeantimerPatternReco::theAlgoName [private]

Definition at line 103 of file DTMeantimerPatternReco.h.

Referenced by algoName().

double DTMeantimerPatternReco::theAlphaMaxPhi [private]

Definition at line 106 of file DTMeantimerPatternReco.h.

Referenced by buildSegments(), and DTMeantimerPatternReco().

double DTMeantimerPatternReco::theAlphaMaxTheta [private]

Definition at line 105 of file DTMeantimerPatternReco.h.

Referenced by buildSegments(), and DTMeantimerPatternReco().

DTSegmentCleaner* DTMeantimerPatternReco::theCleaner [private]

Definition at line 112 of file DTMeantimerPatternReco.h.

Referenced by buildSegments(), DTMeantimerPatternReco(), and ~DTMeantimerPatternReco().

edm::ESHandle<DTGeometry> DTMeantimerPatternReco::theDTGeometry [private]

Definition at line 116 of file DTMeantimerPatternReco.h.

Referenced by initHits(), and setES().

unsigned int DTMeantimerPatternReco::theMaxAllowedHits [private]

Definition at line 104 of file DTMeantimerPatternReco.h.

Referenced by buildSegments(), and DTMeantimerPatternReco().

double DTMeantimerPatternReco::theMaxChi2 [private]

Definition at line 107 of file DTMeantimerPatternReco.h.

Referenced by DTMeantimerPatternReco().

double DTMeantimerPatternReco::theMaxT0 [private]

Definition at line 108 of file DTMeantimerPatternReco.h.

Referenced by DTMeantimerPatternReco().

double DTMeantimerPatternReco::theMinT0 [private]

Definition at line 109 of file DTMeantimerPatternReco.h.

Referenced by DTMeantimerPatternReco().

DTSegmentUpdator* DTMeantimerPatternReco::theUpdator [private]

Definition at line 111 of file DTMeantimerPatternReco.h.

Referenced by DTMeantimerPatternReco(), reconstruct(), setES(), and ~DTMeantimerPatternReco().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:18:57 2009 for CMSSW by  doxygen 1.5.4