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, edm::ConsumesCollector cc)
 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
 
const edm::ESGetToken< DTGeometry, MuonGeometryRecordtheDTGeometryToken
 
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 42 of file DTRefitAndCombineReco4D.h.

Constructor & Destructor Documentation

◆ DTRefitAndCombineReco4D()

DTRefitAndCombineReco4D::DTRefitAndCombineReco4D ( const edm::ParameterSet pset,
edm::ConsumesCollector  cc 
)

Constructor.

Definition at line 30 of file DTRefitAndCombineReco4D.cc.

References gpuPixelDoublets::cc, debug, muonDTDigis_cfi::pset, theMaxChi2forPhi, and theUpdator.

31  : DTRecSegment4DBaseAlgo(pset), theAlgoName("DTRefitAndCombineReco4D"), theDTGeometryToken(cc.esConsumes()) {
32  // debug parameter
33  debug = pset.getUntrackedParameter<bool>("debug");
34 
35  // the updator
37 
38  // the max allowd chi^2 for the fit of th combination of two phi segments
39  theMaxChi2forPhi = pset.getParameter<double>("MaxChi2forPhi");
40 }
DTRecSegment4DBaseAlgo(const edm::ParameterSet &)
Constructor.
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > theDTGeometryToken

◆ ~DTRefitAndCombineReco4D()

DTRefitAndCombineReco4D::~DTRefitAndCombineReco4D ( )
inlineoverride

Destructor.

Definition at line 48 of file DTRefitAndCombineReco4D.h.

48 {}

Member Function Documentation

◆ algoName()

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

Implements DTRecSegment4DBaseAlgo.

Definition at line 53 of file DTRefitAndCombineReco4D.h.

References theAlgoName.

53 { return theAlgoName; }

◆ reconstruct()

OwnVector< DTRecSegment4D > DTRefitAndCombineReco4D::reconstruct ( )
overridevirtual

Operations.

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

Implements DTRecSegment4DBaseAlgo.

Definition at line 83 of file DTRefitAndCombineReco4D.cc.

References gather_cfg::cout, debug, DTChamber::id(), phi, refitSuperSegments(), mps_fire::result, DTChamber::superLayer(), theChamber, theSegments2DTheta, theUpdator, GeomDet::toGlobal(), GeomDet::toLocal(), and DTSegmentUpdator::update().

83  {
85 
86  if (debug)
87  cout << "Segments in " << theChamber->id() << endl;
88 
89  vector<DTChamberRecSegment2D> resultPhi = refitSuperSegments();
90 
91  if (debug)
92  cout << "There are " << resultPhi.size() << " Phi cand" << endl;
93 
94  bool hasZed = false;
95 
96  // has this chamber the Z-superlayer?
97  if (!theSegments2DTheta.empty()) {
98  hasZed = !theSegments2DTheta.empty();
99  if (debug)
100  cout << "There are " << theSegments2DTheta.size() << " Theta cand" << endl;
101  } else {
102  if (debug)
103  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(); phi != resultPhi.end(); ++phi) {
109  if (hasZed) {
110  // Create all the 4D-segment combining the Z view with the Phi one
111  // loop over the Z segments
112  for (vector<DTSLRecSegment2D>::const_iterator zed = theSegments2DTheta.begin(); zed != theSegments2DTheta.end();
113  ++zed) {
114  //>> Important!!
115  DTSuperLayerId ZedSegSLId(zed->geographicalId().rawId());
116 
117  const LocalPoint posZInCh =
118  theChamber->toLocal(theChamber->superLayer(ZedSegSLId)->toGlobal(zed->localPosition()));
119  const LocalVector dirZInCh =
120  theChamber->toLocal(theChamber->superLayer(ZedSegSLId)->toGlobal(zed->localDirection()));
121 
122  DTRecSegment4D* newSeg = new DTRecSegment4D(*phi, *zed, posZInCh, dirZInCh);
123  //<<
124 
126  theUpdator->update(newSeg, false, true);
127  if (debug)
128  cout << "Created a 4D seg " << endl;
129  result.push_back(newSeg);
130  }
131  } else {
132  // Only phi
133  DTRecSegment4D* newSeg = new DTRecSegment4D(*phi);
134  if (debug)
135  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(); zed != theSegments2DTheta.end();
143  ++zed) {
144  // Important!!
145  DTSuperLayerId ZedSegSLId(zed->geographicalId().rawId());
146 
147  const LocalPoint posZInCh =
148  theChamber->toLocal(theChamber->superLayer(ZedSegSLId)->toGlobal(zed->localPosition()));
149  const LocalVector dirZInCh =
150  theChamber->toLocal(theChamber->superLayer(ZedSegSLId)->toGlobal(zed->localDirection()));
151 
152  DTRecSegment4D* newSeg = new DTRecSegment4D(*zed, posZInCh, dirZInCh);
153  //<<
154 
155  if (debug)
156  cout << "Created a 4D segment using only the 2D Theta segment" << endl;
157  result.push_back(newSeg);
158  }
159  }
160  }
161 
162  return result;
163 }
std::vector< DTSLRecSegment2D > theSegments2DTheta
DTChamberId id() const
Return the DTChamberId of this chamber.
Definition: DTChamber.cc:32
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
void update(DTRecSegment4D *seg, const bool calcT0, bool allow3par) const
recompute hits position and refit the segment4D
const DTSuperLayer * superLayer(const DTSuperLayerId &id) const
Return the superlayer corresponding to the given id.
Definition: DTChamber.cc:53
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
std::vector< DTChamberRecSegment2D > refitSuperSegments()

◆ refitSuperSegments()

vector< DTChamberRecSegment2D > DTRefitAndCombineReco4D::refitSuperSegments ( )
private

Definition at line 165 of file DTRefitAndCombineReco4D.cc.

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

Referenced by reconstruct().

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

◆ setChamber()

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:90
edm::ESHandle< DTGeometry > theDTGeometry

◆ setDTRecHit1DContainer()

void DTRefitAndCombineReco4D::setDTRecHit1DContainer ( edm::Handle< DTRecHitCollection all1DHits)
inlineoverridevirtual

Implements DTRecSegment4DBaseAlgo.

Definition at line 57 of file DTRefitAndCombineReco4D.h.

57 {}

◆ setDTRecSegment2DContainer()

void DTRefitAndCombineReco4D::setDTRecSegment2DContainer ( edm::Handle< DTRecSegment2DCollection all2DSegments)
overridevirtual

Implements DTRecSegment4DBaseAlgo.

Definition at line 53 of file DTRefitAndCombineReco4D.cc.

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

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

◆ setES()

void DTRefitAndCombineReco4D::setES ( const edm::EventSetup setup)
overridevirtual

Implements DTRecSegment4DBaseAlgo.

Definition at line 42 of file DTRefitAndCombineReco4D.cc.

References DTSegmentUpdator::setES(), singleTopDQM_cfi::setup, theDTGeometry, theDTGeometryToken, and theUpdator.

42  {
45  // the2DAlgo->setES(setup);
46 }
void setES(const edm::EventSetup &setup)
set the setup
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > theDTGeometryToken
edm::ESHandle< DTGeometry > theDTGeometry

◆ wants2DSegments()

bool DTRefitAndCombineReco4D::wants2DSegments ( )
inlineoverridevirtual

Implements DTRecSegment4DBaseAlgo.

Definition at line 60 of file DTRefitAndCombineReco4D.h.

60 { return true; }

Member Data Documentation

◆ debug

bool DTRefitAndCombineReco4D::debug
private

◆ theAlgoName

std::string DTRefitAndCombineReco4D::theAlgoName
private

Definition at line 66 of file DTRefitAndCombineReco4D.h.

Referenced by algoName().

◆ theChamber

const DTChamber* DTRefitAndCombineReco4D::theChamber
private

Definition at line 83 of file DTRefitAndCombineReco4D.h.

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

◆ theDTGeometry

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

Definition at line 74 of file DTRefitAndCombineReco4D.h.

Referenced by setChamber(), and setES().

◆ theDTGeometryToken

const edm::ESGetToken<DTGeometry, MuonGeometryRecord> DTRefitAndCombineReco4D::theDTGeometryToken
private

Definition at line 75 of file DTRefitAndCombineReco4D.h.

Referenced by setES().

◆ theMaxChi2forPhi

double DTRefitAndCombineReco4D::theMaxChi2forPhi
private

Definition at line 68 of file DTRefitAndCombineReco4D.h.

Referenced by DTRefitAndCombineReco4D(), and refitSuperSegments().

◆ theSegments2DPhi1

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

Definition at line 84 of file DTRefitAndCombineReco4D.h.

Referenced by refitSuperSegments(), and setDTRecSegment2DContainer().

◆ theSegments2DPhi2

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

Definition at line 86 of file DTRefitAndCombineReco4D.h.

Referenced by refitSuperSegments(), and setDTRecSegment2DContainer().

◆ theSegments2DTheta

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

Definition at line 85 of file DTRefitAndCombineReco4D.h.

Referenced by reconstruct(), and setDTRecSegment2DContainer().

◆ theUpdator

DTSegmentUpdator* DTRefitAndCombineReco4D::theUpdator
private