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.

30  : DTRecSegment4DBaseAlgo(pset), theAlgoName("DTRefitAndCombineReco4D") {
31  // debug parameter
32  debug = pset.getUntrackedParameter<bool>("debug");
33 
34  // the updator
35  theUpdator = new DTSegmentUpdator(pset);
36 
37  // the max allowd chi^2 for the fit of th combination of two phi segments
38  theMaxChi2forPhi = pset.getParameter<double>("MaxChi2forPhi");
39 }
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 45 of file DTRefitAndCombineReco4D.h.

45 {};

Member Function Documentation

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

Implements DTRecSegment4DBaseAlgo.

Definition at line 50 of file DTRefitAndCombineReco4D.h.

References singleTopDQM_cfi::setup.

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

Operations.

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

Implements DTRecSegment4DBaseAlgo.

Definition at line 82 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().

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

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

Referenced by reconstruct().

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

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

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

Implements DTRecSegment4DBaseAlgo.

Definition at line 54 of file DTRefitAndCombineReco4D.h.

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

Implements DTRecSegment4DBaseAlgo.

Definition at line 52 of file DTRefitAndCombineReco4D.cc.

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

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

Implements DTRecSegment4DBaseAlgo.

Definition at line 41 of file DTRefitAndCombineReco4D.cc.

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

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

Implements DTRecSegment4DBaseAlgo.

Definition at line 57 of file DTRefitAndCombineReco4D.h.

57 { return true; }

Member Data Documentation

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

Definition at line 63 of file DTRefitAndCombineReco4D.h.

const DTChamber* DTRefitAndCombineReco4D::theChamber
private

Definition at line 79 of file DTRefitAndCombineReco4D.h.

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

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

Definition at line 71 of file DTRefitAndCombineReco4D.h.

Referenced by setChamber(), and setES().

double DTRefitAndCombineReco4D::theMaxChi2forPhi
private

Definition at line 65 of file DTRefitAndCombineReco4D.h.

Referenced by DTRefitAndCombineReco4D(), and refitSuperSegments().

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

Definition at line 80 of file DTRefitAndCombineReco4D.h.

Referenced by refitSuperSegments(), and setDTRecSegment2DContainer().

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

Definition at line 82 of file DTRefitAndCombineReco4D.h.

Referenced by refitSuperSegments(), and setDTRecSegment2DContainer().

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

Definition at line 81 of file DTRefitAndCombineReco4D.h.

Referenced by reconstruct(), and setDTRecSegment2DContainer().

DTSegmentUpdator* DTRefitAndCombineReco4D::theUpdator
private