CMS 3D CMS Logo

DTSegment2DQuality.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author S. Bolognesi and G. Cerminara - INFN Torino
5  */
6 
7 #include <iostream>
8 #include <map>
9 
20 
21 #include "DTSegment2DQuality.h"
22 #include "Histograms.h"
23 
24 using namespace std;
25 using namespace edm;
26 
27 namespace dtsegment2d {
28  struct Histograms {
29  std::unique_ptr<HRes2DHit> h2DHitRPhi;
30  std::unique_ptr<HRes2DHit> h2DHitRZ;
31  std::unique_ptr<HRes2DHit> h2DHitRZ_W0;
32  std::unique_ptr<HRes2DHit> h2DHitRZ_W1;
33  std::unique_ptr<HRes2DHit> h2DHitRZ_W2;
34 
35  std::unique_ptr<HEff2DHit> h2DHitEff_RPhi;
36  std::unique_ptr<HEff2DHit> h2DHitEff_RZ;
37  std::unique_ptr<HEff2DHit> h2DHitEff_RZ_W0;
38  std::unique_ptr<HEff2DHit> h2DHitEff_RZ_W1;
39  std::unique_ptr<HEff2DHit> h2DHitEff_RZ_W2;
40  };
41 } // namespace dtsegment2d
42 
43 using namespace dtsegment2d;
44 
45 // Constructor
47  // get the debug parameter for verbose output
48  debug_ = pset.getUntrackedParameter<bool>("debug");
50  // the name of the simhit collection
51  simHitLabel_ = pset.getUntrackedParameter<InputTag>("simHitLabel");
52  simHitToken_ = consumes<PSimHitContainer>(pset.getUntrackedParameter<InputTag>("simHitLabel"));
53  // the name of the 2D rec hit collection
54  segment2DLabel_ = pset.getUntrackedParameter<InputTag>("segment2DLabel");
55  segment2DToken_ = consumes<DTRecSegment2DCollection>(pset.getUntrackedParameter<InputTag>("segment2DLabel"));
56 
57  // sigma resolution on position
58  sigmaResPos_ = pset.getParameter<double>("sigmaResPos");
59  // sigma resolution on angle
60  sigmaResAngle_ = pset.getParameter<double>("sigmaResAngle");
61 
62  if (debug_) {
63  cout << "[DTSegment2DQuality] Constructor called " << endl;
64  }
65 }
66 
68  edm::Run const &run,
69  edm::EventSetup const &setup,
70  Histograms &histograms) const {
71  histograms.h2DHitRPhi = std::make_unique<HRes2DHit>("RPhi", booker, true, true);
72  histograms.h2DHitRZ = std::make_unique<HRes2DHit>("RZ", booker, true, true);
73  histograms.h2DHitRZ_W0 = std::make_unique<HRes2DHit>("RZ_W0", booker, true, true);
74  histograms.h2DHitRZ_W1 = std::make_unique<HRes2DHit>("RZ_W1", booker, true, true);
75  histograms.h2DHitRZ_W2 = std::make_unique<HRes2DHit>("RZ_W2", booker, true, true);
76 
77  histograms.h2DHitEff_RPhi = std::make_unique<HEff2DHit>("RPhi", booker);
78  histograms.h2DHitEff_RZ = std::make_unique<HEff2DHit>("RZ", booker);
79  histograms.h2DHitEff_RZ_W0 = std::make_unique<HEff2DHit>("RZ_W0", booker);
80  histograms.h2DHitEff_RZ_W1 = std::make_unique<HEff2DHit>("RZ_W1", booker);
81  histograms.h2DHitEff_RZ_W2 = std::make_unique<HEff2DHit>("RZ_W2", booker);
82  if (debug_) {
83  cout << "[DTSegment2DQuality] hitsos created " << endl;
84  }
85 }
86 
87 // The real analysis
89  edm::EventSetup const &setup,
90  Histograms const &histograms) const {
91  // Get the DT Geometry
92  const DTGeometry &dtGeom = setup.getData(muonGeomToken_);
93 
94  // Get the SimHit collection from the event
96  event.getByToken(simHitToken_, simHits); // FIXME: second string to be removed
97 
98  // Map simHits by sl
99  map<DTSuperLayerId, PSimHitContainer> simHitsPerSl;
100  for (const auto &simHit : *simHits) {
101  // Create the id of the sl (the simHits in the DT known their wireId)
102  DTSuperLayerId slId = ((DTWireId(simHit.detUnitId())).layerId()).superlayerId();
103  // Fill the map
104  simHitsPerSl[slId].push_back(simHit);
105  }
106 
107  // Get the 2D rechits from the event
109  event.getByToken(segment2DToken_, segment2Ds);
110 
111  if (not segment2Ds.isValid()) {
112  if (debug_) {
113  cout << "[DTSegment2DQuality]**Warning: no 2DSegments with label: " << segment2DLabel_
114  << " in this event, skipping !" << endl;
115  }
116  return;
117  }
118 
119  // Loop over all superlayers containing a segment
121  for (slId = segment2Ds->id_begin(); slId != segment2Ds->id_end(); ++slId) {
122  //------------------------- simHits ---------------------------//
123  // Get simHits of each superlayer
124  PSimHitContainer simHits = simHitsPerSl[(*slId)];
125 
126  // Map simhits per wire
127  map<DTWireId, PSimHitContainer> simHitsPerWire = DTHitQualityUtils::mapSimHitsPerWire(simHits);
128  map<DTWireId, const PSimHit *> muSimHitPerWire = DTHitQualityUtils::mapMuSimHitsPerWire(simHitsPerWire);
129  int nMuSimHit = muSimHitPerWire.size();
130  if (nMuSimHit == 0 or nMuSimHit == 1) {
131  if (debug_ and nMuSimHit == 1) {
132  cout << "[DTSegment2DQuality] Only " << nMuSimHit << " mu SimHit in this SL, skipping !" << endl;
133  }
134  continue; // If no or only one mu SimHit is found skip this SL
135  }
136  if (debug_) {
137  cout << "=== SL " << (*slId) << " has " << nMuSimHit << " SimHits" << endl;
138  }
139 
140  // Find outer and inner mu SimHit to build a segment
141  pair<const PSimHit *, const PSimHit *> inAndOutSimHit = DTHitQualityUtils::findMuSimSegment(muSimHitPerWire);
142  // Check that outermost and innermost SimHit are not the same
143  if (inAndOutSimHit.first == inAndOutSimHit.second) {
144  cout << "[DTHitQualityUtils]***Warning: outermost and innermost SimHit "
145  "are the same !"
146  << endl;
147  continue;
148  }
149 
150  // Find direction and position of the sim Segment in SL RF
151  pair<LocalVector, LocalPoint> dirAndPosSimSegm =
152  DTHitQualityUtils::findMuSimSegmentDirAndPos(inAndOutSimHit, (*slId), dtGeom);
153 
154  LocalVector simSegmLocalDir = dirAndPosSimSegm.first;
155  LocalPoint simSegmLocalPos = dirAndPosSimSegm.second;
156  if (debug_) {
157  cout << " Simulated segment: local direction " << simSegmLocalDir << endl
158  << " local position " << simSegmLocalPos << endl;
159  }
160  const DTSuperLayer *superLayer = dtGeom.superLayer(*slId);
161  GlobalPoint simSegmGlobalPos = superLayer->toGlobal(simSegmLocalPos);
162 
163  // Atan(x/z) angle and x position in SL RF
164  float angleSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).first;
165  float posSimSeg = simSegmLocalPos.x();
166  // Position (in eta, phi coordinates) in the global RF
167  float etaSimSeg = simSegmGlobalPos.eta();
168  float phiSimSeg = simSegmGlobalPos.phi();
169 
170  //---------------------------- recHits --------------------------//
171  // Get the range of rechit for the corresponding slId
172  bool recHitFound = false;
173  DTRecSegment2DCollection::range range = segment2Ds->get(*slId);
174  int nsegm = distance(range.first, range.second);
175  if (debug_) {
176  cout << " Sl: " << *slId << " has " << nsegm << " 2D segments" << endl;
177  }
178 
179  if (nsegm != 0) {
180  // Find the best RecHit: look for the 2D RecHit with the angle closest
181  // to that of segment made of SimHits.
182  // RecHits must have delta alpha and delta position within 5 sigma of
183  // the residual distribution (we are looking for residuals of segments
184  // usefull to the track fit) for efficency purpose
185  const DTRecSegment2D *bestRecHit = nullptr;
186  bool bestRecHitFound = false;
187  double deltaAlpha = 99999;
188 
189  // Loop over the recHits of this slId
190  for (DTRecSegment2DCollection::const_iterator segment2D = range.first; segment2D != range.second; ++segment2D) {
191  // Check the dimension
192  if ((*segment2D).dimension() != 2) {
193  if (debug_) {
194  cout << "[DTSegment2DQuality]***Error: This is not 2D segment !!!" << endl;
195  }
196  abort();
197  }
198  // Segment Local Direction and position (in SL RF)
199  LocalVector recSegDirection = (*segment2D).localDirection();
200  LocalPoint recSegPosition = (*segment2D).localPosition();
201 
202  float recSegAlpha = DTHitQualityUtils::findSegmentAlphaAndBeta(recSegDirection).first;
203  if (debug_) {
204  cout << " RecSegment direction: " << recSegDirection << endl
205  << " position : " << recSegPosition << endl
206  << " alpha : " << recSegAlpha << endl;
207  }
208 
209  if (fabs(recSegAlpha - angleSimSeg) < deltaAlpha) {
210  deltaAlpha = fabs(recSegAlpha - angleSimSeg);
211  bestRecHit = &(*segment2D);
212  bestRecHitFound = true;
213  }
214  } // End of Loop over all 2D RecHits
215 
216  if (bestRecHitFound) {
217  // Best rechit direction and position in SL RF
218  LocalPoint bestRecHitLocalPos = bestRecHit->localPosition();
219  LocalVector bestRecHitLocalDir = bestRecHit->localDirection();
220 
221  LocalError bestRecHitLocalPosErr = bestRecHit->localPositionError();
222  LocalError bestRecHitLocalDirErr = bestRecHit->localDirectionError();
223 
224  float angleBestRHit = DTHitQualityUtils::findSegmentAlphaAndBeta(bestRecHitLocalDir).first;
225 
226  if (fabs(angleBestRHit - angleSimSeg) < 5 * sigmaResAngle_ and
227  fabs(bestRecHitLocalPos.x() - posSimSeg) < 5 * sigmaResPos_) {
228  recHitFound = true;
229  }
230 
231  // Fill Residual histos
232  HRes2DHit *hRes = nullptr;
233  if ((*slId).superlayer() == 1 or (*slId).superlayer() == 3) { // RPhi SL
234  hRes = histograms.h2DHitRPhi.get();
235  } else if ((*slId).superlayer() == 2) { // RZ SL
236  histograms.h2DHitRZ->fill(angleSimSeg,
237  angleBestRHit,
238  posSimSeg,
239  bestRecHitLocalPos.x(),
240  etaSimSeg,
241  phiSimSeg,
242  sqrt(bestRecHitLocalPosErr.xx()),
243  sqrt(bestRecHitLocalDirErr.xx()));
244  if (abs((*slId).wheel()) == 0) {
245  hRes = histograms.h2DHitRZ_W0.get();
246  } else if (abs((*slId).wheel()) == 1) {
247  hRes = histograms.h2DHitRZ_W1.get();
248  } else if (abs((*slId).wheel()) == 2) {
249  hRes = histograms.h2DHitRZ_W2.get();
250  }
251  }
252  hRes->fill(angleSimSeg,
253  angleBestRHit,
254  posSimSeg,
255  bestRecHitLocalPos.x(),
256  etaSimSeg,
257  phiSimSeg,
258  sqrt(bestRecHitLocalPosErr.xx()),
259  sqrt(bestRecHitLocalDirErr.xx()));
260  }
261  } // end of if (nsegm != 0)
262 
263  // Fill Efficiency plot
264  HEff2DHit *hEff = nullptr;
265  if ((*slId).superlayer() == 1 or (*slId).superlayer() == 3) { // RPhi SL
266  hEff = histograms.h2DHitEff_RPhi.get();
267  } else if ((*slId).superlayer() == 2) { // RZ SL
268  histograms.h2DHitEff_RZ->fill(etaSimSeg, phiSimSeg, posSimSeg, angleSimSeg, recHitFound);
269  if (abs((*slId).wheel()) == 0) {
270  hEff = histograms.h2DHitEff_RZ_W0.get();
271  } else if (abs((*slId).wheel()) == 1) {
272  hEff = histograms.h2DHitEff_RZ_W1.get();
273  } else if (abs((*slId).wheel()) == 2) {
274  hEff = histograms.h2DHitEff_RZ_W2.get();
275  }
276  }
277  hEff->fill(etaSimSeg, phiSimSeg, posSimSeg, angleSimSeg, recHitFound);
278  } // End of loop over superlayers
279 }
280 
281 // declare this as a framework plugin
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &, dtsegment2d::Histograms &) const override
Book the DQM plots.
std::pair< double, double > findSegmentAlphaAndBeta(const LocalVector &direction)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
std::map< DTWireId, edm::PSimHitContainer > mapSimHitsPerWire(const edm::PSimHitContainer &simhits)
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
LocalPoint localPosition() const override
local position in SL frame
edm::EDGetTokenT< edm::PSimHitContainer > simHitToken_
LocalVector localDirection() const override
the local direction in SL frame
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
void fill(float etaSimSegm, float phiSimSegm, float posSimSegm, float angleSimSegm, bool fillRecHit)
Definition: Histograms.h:354
T eta() const
Definition: PV3DBase.h:73
std::pair< LocalVector, LocalPoint > findMuSimSegmentDirAndPos(const std::pair< const PSimHit *, const PSimHit *> &inAndOutSimHit, const DetId detId, const DTGeometry &muonGeom)
identifier iterator
Definition: RangeMap.h:130
edm::ESGetToken< DTGeometry, MuonGeometryRecord > muonGeomToken_
DTSegment2DQuality(const edm::ParameterSet &pset)
Constructor.
std::atomic< bool > debug
std::unique_ptr< HRes2DHit > h2DHitRZ_W2
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
T x() const
Definition: PV3DBase.h:59
const DTSuperLayer * superLayer(const DTSuperLayerId &id) const
Return a DTSuperLayer given its id.
Definition: DTGeometry.cc:92
T sqrt(T t)
Definition: SSEVec.h:23
LocalError localDirectionError() const override
the local direction error (xx,xy,yy) in SL frame: only xx is not 0.
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
std::map< DTWireId, const PSimHit * > mapMuSimHitsPerWire(const std::map< DTWireId, edm::PSimHitContainer > &simHitWireMap)
Create a map between the Mu SimHits and corresponding MuBarWireId ;.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
LocalError localPositionError() const override
local position error in SL frame
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< DTRecSegment2DCollection > segment2DToken_
std::unique_ptr< HEff2DHit > h2DHitEff_RZ_W1
std::unique_ptr< HRes2DHit > h2DHitRZ_W0
std::unique_ptr< HEff2DHit > h2DHitEff_RZ_W0
edm::InputTag simHitLabel_
std::unique_ptr< HEff2DHit > h2DHitEff_RPhi
std::unique_ptr< HRes2DHit > h2DHitRPhi
std::pair< const PSimHit *, const PSimHit * > findMuSimSegment(const std::map< DTWireId, const PSimHit *> &mapWireAndMuSimHit)
bool isValid() const
Definition: HandleBase.h:70
void dqmAnalyze(edm::Event const &, edm::EventSetup const &, dtsegment2d::Histograms const &) const override
Perform the real analysis.
HLT enums.
std::unique_ptr< HRes2DHit > h2DHitRZ_W1
std::vector< PSimHit > PSimHitContainer
std::unique_ptr< HEff2DHit > h2DHitEff_RZ
edm::InputTag segment2DLabel_
float xx() const
Definition: LocalError.h:22
Definition: event.py:1
Definition: Run.h:45
std::unique_ptr< HEff2DHit > h2DHitEff_RZ_W2
std::unique_ptr< HRes2DHit > h2DHitRZ