CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
DTRefitAndCombineReco4D Class Reference

#include <DTRefitAndCombineReco4D.h>

Inheritance diagram for DTRefitAndCombineReco4D:
DTRecSegment4DBaseAlgo

Public Member Functions

std::string algoName () const override
 
 DTRefitAndCombineReco4D (const edm::ParameterSet &pset)
 Constructor. More...
 
edm::OwnVector< DTRecSegment4Dreconstruct () override
 Operations. More...
 
void setChamber (const DTChamberId &chId) override
 
void setDTRecHit1DContainer (edm::Handle< DTRecHitCollection > all1DHits) override
 
void setDTRecSegment2DContainer (edm::Handle< DTRecSegment2DCollection > all2DSegments) override
 
void setES (const edm::EventSetup &setup) override
 
bool wants2DSegments () override
 
 ~DTRefitAndCombineReco4D () override
 Destructor. More...
 
- Public Member Functions inherited from DTRecSegment4DBaseAlgo
 DTRecSegment4DBaseAlgo (const edm::ParameterSet &)
 Constructor. More...
 
virtual ~DTRecSegment4DBaseAlgo ()
 Destructor. More...
 

Private Member Functions

std::vector< DTChamberRecSegment2DrefitSuperSegments ()
 

Private Attributes

bool debug
 
std::string theAlgoName
 
const DTChambertheChamber
 
edm::ESHandle< DTGeometrytheDTGeometry
 
double theMaxChi2forPhi
 
std::vector< DTSLRecSegment2DtheSegments2DPhi1
 
std::vector< DTSLRecSegment2DtheSegments2DPhi2
 
std::vector< DTSLRecSegment2DtheSegments2DTheta
 
DTSegmentUpdatortheUpdator
 

Detailed Description

Algo for reconstructing 4d segment in DT refitting the 2D phi SL hits and combining the results with the theta view.

Author
Stefano Lacaprara - INFN Legnaro stefa.nosp@m.no.l.nosp@m.acapr.nosp@m.ara@.nosp@m.pd.in.nosp@m.fn.i.nosp@m.t
Riccardo Bellan - INFN TO ricca.nosp@m.rdo..nosp@m.bella.nosp@m.n@ce.nosp@m.rn.ch

Definition at line 39 of file DTRefitAndCombineReco4D.h.

Constructor & Destructor Documentation

DTRefitAndCombineReco4D::DTRefitAndCombineReco4D ( const edm::ParameterSet pset)

Constructor.

Definition at line 29 of file DTRefitAndCombineReco4D.cc.

References debug, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), theMaxChi2forPhi, and theUpdator.

29  :
30 DTRecSegment4DBaseAlgo(pset), theAlgoName("DTRefitAndCombineReco4D"){
31 
32  // debug parameter
33  debug = pset.getUntrackedParameter<bool>("debug");
34 
35  // the updator
36  theUpdator = new DTSegmentUpdator(pset);
37 
38  // the max allowd chi^2 for the fit of th combination of two phi segments
39  theMaxChi2forPhi = pset.getParameter<double>("MaxChi2forPhi");
40 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
DTRecSegment4DBaseAlgo(const edm::ParameterSet &)
Constructor.
DTRefitAndCombineReco4D::~DTRefitAndCombineReco4D ( )
inlineoverride

Destructor.

Definition at line 47 of file DTRefitAndCombineReco4D.h.

47 {};

Member Function Documentation

std::string DTRefitAndCombineReco4D::algoName ( ) const
inlineoverridevirtual

Implements DTRecSegment4DBaseAlgo.

Definition at line 53 of file DTRefitAndCombineReco4D.h.

References GeneralSetup::setup().

53 { return theAlgoName; }
OwnVector< DTRecSegment4D > DTRefitAndCombineReco4D::reconstruct ( )
overridevirtual

Operations.

4d segment: I have the pos along the wire => further update!

Implements DTRecSegment4DBaseAlgo.

Definition at line 87 of file DTRefitAndCombineReco4D.cc.

References gather_cfg::cout, debug, DTChamber::id(), phi, edm::OwnVector< T, P >::push_back(), refitSuperSegments(), mps_fire::result, DTChamber::superLayer(), theChamber, theSegments2DTheta, theUpdator, GeomDet::toGlobal(), GeomDet::toLocal(), and DTSegmentUpdator::update().

87  {
89 
90  if (debug) cout << "Segments in " << theChamber->id() << endl;
91 
92  vector<DTChamberRecSegment2D> resultPhi = refitSuperSegments();
93 
94  if (debug) cout << "There are " << resultPhi.size() << " Phi cand" << endl;
95 
96  bool hasZed=false;
97 
98  // has this chamber the Z-superlayer?
99  if (!theSegments2DTheta.empty()){
100  hasZed = !theSegments2DTheta.empty();
101  if (debug) cout << "There are " << theSegments2DTheta.size() << " Theta cand" << endl;
102  } else {
103  if (debug) cout << "No Theta SL" << endl;
104  }
105 
106  // Now I want to build the concrete DTRecSegment4D.
107  if (!resultPhi.empty()) {
108  for (vector<DTChamberRecSegment2D>::const_iterator phi=resultPhi.begin();
109  phi!=resultPhi.end(); ++phi) {
110 
111  if (hasZed) {
112 
113  // Create all the 4D-segment combining the Z view with the Phi one
114  // loop over the Z segments
115  for(vector<DTSLRecSegment2D>::const_iterator zed = theSegments2DTheta.begin();
116  zed != theSegments2DTheta.end(); ++zed){
117 
118  //>> Important!!
119  DTSuperLayerId ZedSegSLId(zed->geographicalId().rawId());
120 
121  const LocalPoint posZInCh = theChamber->toLocal( theChamber->superLayer(ZedSegSLId)->toGlobal(zed->localPosition() )) ;
122  const LocalVector dirZInCh = theChamber->toLocal( theChamber->superLayer(ZedSegSLId)->toGlobal(zed->localDirection() )) ;
123 
124  DTRecSegment4D* newSeg = new DTRecSegment4D(*phi,*zed,posZInCh,dirZInCh);
125  //<<
126 
128  theUpdator->update(newSeg,false,true);
129  if (debug) cout << "Created a 4D seg " << endl;
130  result.push_back(newSeg);
131  }
132  } else {
133  // Only phi
134  DTRecSegment4D* newSeg = new DTRecSegment4D(*phi);
135  if (debug) cout << "Created a 4D segment using only the 2D Phi segment" << endl;
136  result.push_back(newSeg);
137  }
138  }
139  } else {
140  // DTRecSegment4D from zed projection only (unlikely, not so useful, but...)
141  if (hasZed) {
142  for(vector<DTSLRecSegment2D>::const_iterator zed = theSegments2DTheta.begin();
143  zed != theSegments2DTheta.end(); ++zed){
144 
145  // Important!!
146  DTSuperLayerId ZedSegSLId(zed->geographicalId().rawId());
147 
148  const LocalPoint posZInCh = theChamber->toLocal( theChamber->superLayer(ZedSegSLId)->toGlobal(zed->localPosition() )) ;
149  const LocalVector dirZInCh = theChamber->toLocal( theChamber->superLayer(ZedSegSLId)->toGlobal(zed->localDirection() )) ;
150 
151  DTRecSegment4D* newSeg = new DTRecSegment4D( *zed,posZInCh,dirZInCh);
152  //<<
153 
154  if (debug) cout << "Created a 4D segment using only the 2D Theta segment" << endl;
155  result.push_back(newSeg);
156  }
157  }
158  }
159 
160  return result;
161 }
const DTSuperLayer * superLayer(const DTSuperLayerId &id) const
Return the superlayer corresponding to the given id.
Definition: DTChamber.cc:65
std::vector< DTSLRecSegment2D > theSegments2DTheta
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:69
void push_back(D *&d)
Definition: OwnVector.h:290
DTChamberId id() const
Return the DTChamberId of this chamber.
Definition: DTChamber.cc:33
std::vector< DTChamberRecSegment2D > refitSuperSegments()
void update(DTRecSegment4D *seg, const bool calcT0, bool allow3par) const
recompute hits position and refit the segment4D
vector< DTChamberRecSegment2D > DTRefitAndCombineReco4D::refitSuperSegments ( )
private

Definition at line 163 of file DTRefitAndCombineReco4D.cc.

References DTRecSegment2D::chi2(), popcon2dropbox::copy(), DTSegmentUpdator::fit(), mps_fire::result, theMaxChi2forPhi, theSegments2DPhi1, theSegments2DPhi2, and theUpdator.

Referenced by reconstruct().

163  {
164  vector<DTChamberRecSegment2D> result;
165 
166  //double-loop over all the DTSLRecSegment2D in order to make all the possible pairs
167  for(vector<DTSLRecSegment2D>::const_iterator segment2DPhi1 = theSegments2DPhi1.begin();
168  segment2DPhi1 != theSegments2DPhi1.end(); ++segment2DPhi1){
169  for(vector<DTSLRecSegment2D>::const_iterator segment2DPhi2 = theSegments2DPhi2.begin();
170  segment2DPhi2 != theSegments2DPhi2.end(); ++segment2DPhi2){
171 
172  // check the id
173  if(segment2DPhi1->chamberId() != segment2DPhi2->chamberId())
174  throw cms::Exception("refitSuperSegments")
175  <<"he phi segments have different chamber id"<<std::endl;
176 
177  // create a super phi starting from 2 phi
178  vector<DTRecHit1D> recHitsSeg2DPhi1 = segment2DPhi1->specificRecHits();
179  vector<DTRecHit1D> recHitsSeg2DPhi2 = segment2DPhi2->specificRecHits();
180  // copy the recHitsSeg2DPhi2 in the recHitsSeg2DPhi1 container
181  copy(recHitsSeg2DPhi2.begin(),recHitsSeg2DPhi2.end(),back_inserter(recHitsSeg2DPhi1));
182 
183  const DTChamberId chId = segment2DPhi1->chamberId();
184 
185  // create the super phi
186  DTChamberRecSegment2D superPhi(chId,recHitsSeg2DPhi1);
187 
188  // refit it!
189  theUpdator->fit(&superPhi,false,false);
190 
191  // cut on the chi^2
192  if (superPhi.chi2() > theMaxChi2forPhi)
193  result.push_back(superPhi);
194  }
195  }
196  // TODO clean the container!!!
197  // there are some possible repetition!
198  // maybe using the cleaner, previous a conversion from DTChamberRecSegment2D to DTSegmentCandidate
199  return result;
200 }
def copy(args, dbName)
std::vector< DTSLRecSegment2D > theSegments2DPhi1
std::vector< DTSLRecSegment2D > theSegments2DPhi2
bool fit(DTSegmentCand *seg, bool allow3par, const bool fitdebug) const
void DTRefitAndCombineReco4D::setChamber ( const DTChamberId chId)
overridevirtual

Implements DTRecSegment4DBaseAlgo.

Definition at line 48 of file DTRefitAndCombineReco4D.cc.

References DTGeometry::chamber(), theChamber, and theDTGeometry.

48  {
49  // Set the chamber
51 }
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:117
edm::ESHandle< DTGeometry > theDTGeometry
void DTRefitAndCombineReco4D::setDTRecHit1DContainer ( edm::Handle< DTRecHitCollection all1DHits)
inlineoverridevirtual

Implements DTRecSegment4DBaseAlgo.

Definition at line 57 of file DTRefitAndCombineReco4D.h.

57 {};
void DTRefitAndCombineReco4D::setDTRecSegment2DContainer ( edm::Handle< DTRecSegment2DCollection all2DSegments)
overridevirtual

Implements DTRecSegment4DBaseAlgo.

Definition at line 54 of file DTRefitAndCombineReco4D.cc.

References gather_cfg::cout, debug, DTChamber::id(), theChamber, theSegments2DPhi1, theSegments2DPhi2, and theSegments2DTheta.

54  {
55  theSegments2DPhi1.clear();
56  theSegments2DTheta.clear();
57  theSegments2DPhi2.clear();
58 
59  // Get the chamber
60  // const DTChamber *chamber = theDTGeometry->chamber(chId);
61 
62  const DTChamberId chId = theChamber->id();
63 
64  //Extract the DTRecSegment2DCollection ranges for the three different SL
65  DTRecSegment2DCollection::range rangePhi1 = allHits->get(DTSuperLayerId(chId,1));
66  DTRecSegment2DCollection::range rangeTheta = allHits->get(DTSuperLayerId(chId,2));
67  DTRecSegment2DCollection::range rangePhi2 = allHits->get(DTSuperLayerId(chId,3));
68 
69  // Fill the DTSLRecSegment2D containers for the three different SL
70  vector<DTSLRecSegment2D> segments2DPhi1(rangePhi1.first,rangePhi1.second);
71  vector<DTSLRecSegment2D> segments2DTheta(rangeTheta.first,rangeTheta.second);
72  vector<DTSLRecSegment2D> segments2DPhi2(rangePhi2.first,rangePhi2.second);
73 
74  if(debug)
75  cout << "Number of 2D-segments in the first SL (Phi)" << segments2DPhi1.size() << endl
76  << "Number of 2D-segments in the second SL (Theta)" << segments2DTheta.size() << endl
77  << "Number of 2D-segments in the third SL (Phi)" << segments2DPhi2.size() << endl;
78 
79  theSegments2DPhi1 = segments2DPhi1;
80  theSegments2DTheta = segments2DTheta;
81  theSegments2DPhi2 = segments2DPhi2;
82 }
std::vector< DTSLRecSegment2D > theSegments2DTheta
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
DTChamberId id() const
Return the DTChamberId of this chamber.
Definition: DTChamber.cc:33
std::vector< DTSLRecSegment2D > theSegments2DPhi1
std::vector< DTSLRecSegment2D > theSegments2DPhi2
void DTRefitAndCombineReco4D::setES ( const edm::EventSetup setup)
overridevirtual

Implements DTRecSegment4DBaseAlgo.

Definition at line 42 of file DTRefitAndCombineReco4D.cc.

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

42  {
43  setup.get<MuonGeometryRecord>().get(theDTGeometry);
44  theUpdator->setES(setup);
45  // the2DAlgo->setES(setup);
46 }
void setES(const edm::EventSetup &setup)
set the setup
T get() const
Definition: EventSetup.h:68
edm::ESHandle< DTGeometry > theDTGeometry
bool DTRefitAndCombineReco4D::wants2DSegments ( )
inlineoverridevirtual

Implements DTRecSegment4DBaseAlgo.

Definition at line 60 of file DTRefitAndCombineReco4D.h.

60 {return true;}

Member Data Documentation

bool DTRefitAndCombineReco4D::debug
private
std::string DTRefitAndCombineReco4D::theAlgoName
private

Definition at line 67 of file DTRefitAndCombineReco4D.h.

const DTChamber* DTRefitAndCombineReco4D::theChamber
private

Definition at line 83 of file DTRefitAndCombineReco4D.h.

Referenced by reconstruct(), setChamber(), and setDTRecSegment2DContainer().

edm::ESHandle<DTGeometry> DTRefitAndCombineReco4D::theDTGeometry
private

Definition at line 75 of file DTRefitAndCombineReco4D.h.

Referenced by setChamber(), and setES().

double DTRefitAndCombineReco4D::theMaxChi2forPhi
private

Definition at line 69 of file DTRefitAndCombineReco4D.h.

Referenced by DTRefitAndCombineReco4D(), and refitSuperSegments().

std::vector<DTSLRecSegment2D> DTRefitAndCombineReco4D::theSegments2DPhi1
private

Definition at line 84 of file DTRefitAndCombineReco4D.h.

Referenced by refitSuperSegments(), and setDTRecSegment2DContainer().

std::vector<DTSLRecSegment2D> DTRefitAndCombineReco4D::theSegments2DPhi2
private

Definition at line 86 of file DTRefitAndCombineReco4D.h.

Referenced by refitSuperSegments(), and setDTRecSegment2DContainer().

std::vector<DTSLRecSegment2D> DTRefitAndCombineReco4D::theSegments2DTheta
private

Definition at line 85 of file DTRefitAndCombineReco4D.h.

Referenced by reconstruct(), and setDTRecSegment2DContainer().

DTSegmentUpdator* DTRefitAndCombineReco4D::theUpdator
private