CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
DTChamberEfficiencyTask.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author G. Mila - INFN Torino
5  */
6 
8 
9 //Framework
13 
16 
17 #include <iterator>
18 #include <iostream>
19 #include <cmath>
20 using namespace edm;
21 using namespace std;
22 
24  : muonGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
25  debug = pset.getUntrackedParameter<bool>("debug", false);
26 
27  edm::LogVerbatim("DTDQM|DTMonitorModule|DTChamberEfficiencyTask") << "[DTChamberEfficiencyTask] Constructor called!";
28 
29  parameters = pset;
30 
31  // the name of the 4D rec hits collection
33  consumes<DTRecSegment4DCollection>(edm::InputTag(parameters.getParameter<string>("recHits4DLabel")));
34 
35  // parameters to use for the segment quality check
36  theMinHitsSegment = static_cast<unsigned int>(parameters.getParameter<int>("minHitsSegment"));
37  theMinChi2NormSegment = parameters.getParameter<double>("minChi2NormSegment");
38  // parameter to use for the exstrapolated segment check
39  theMinCloseDist = parameters.getParameter<double>("minCloseDist");
40 
41  // the running modality
42  onlineMonitor = parameters.getUntrackedParameter<bool>("onlineMonitor");
43 
44  // the analysis mode
45  detailedAnalysis = parameters.getUntrackedParameter<bool>("detailedAnalysis");
46 }
47 
49  edm::LogVerbatim("DTDQM|DTMonitorModule|DTChamberEfficiencyTask") << "[DTChamberEfficiencyTask] Destructor called!";
50 }
51 
53  edm::LogVerbatim("DTDQM|DTMonitorModule|DTChamberEfficiencyTask")
54  << "[DTChamberEfficiencyTask]: Begin of LS transition";
55 
56  if (lumiSeg.id().luminosityBlock() % parameters.getUntrackedParameter<int>("ResetCycle", 3) == 0 && onlineMonitor) {
57  for (map<DTChamberId, vector<MonitorElement*> >::const_iterator histo = histosPerCh.begin();
58  histo != histosPerCh.end();
59  histo++) {
60  int size = (*histo).second.size();
61  for (int i = 0; i < size; i++) {
62  (*histo).second[i]->Reset();
63  }
64  }
65  }
66 }
67 
69  // Get the DT Geometry
70  dtGeom = &setup.getData(muonGeomToken_);
71 }
72 
74  edm::Run const& iRun,
75  edm::EventSetup const& context) {
76  ibooker.setCurrentFolder("DT/DTChamberEfficiencyTask");
77 
78  // Loop over all the chambers
79  vector<const DTChamber*>::const_iterator ch_it = dtGeom->chambers().begin();
80  vector<const DTChamber*>::const_iterator ch_end = dtGeom->chambers().end();
81  for (; ch_it != ch_end; ++ch_it) {
82  // histo booking
83  bookHistos(ibooker, (*ch_it)->id());
84  }
85 }
86 
87 // Book a set of histograms for a given Layer
89  edm::LogVerbatim("DTDQM|DTMonitorModule|DTChamberEfficiencyTask") << " Booking histos for CH : " << chId;
90 
91  // Compose the chamber name
92  stringstream wheel;
93  wheel << chId.wheel();
94  stringstream station;
95  station << chId.station();
96  stringstream sector;
97  sector << chId.sector();
98 
99  string HistoName = "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str();
100 
101  ibooker.setCurrentFolder("DT/01-DTChamberEfficiency/Task/Wheel" + wheel.str() + "/Sector" + sector.str() +
102  "/Station" + station.str());
103 
104  // Create the monitor elements
105  vector<MonitorElement*> histos;
106 
107  //efficiency selection cuts
108  // a- number of segments of the top chamber > 0 && number of segments of the bottom chamber > 0
109  // b- number of segments of the middle chamber > 0
110  // c- check of the top and bottom segment quality
111  // d- check if interpolation falls inside the middle chamber
112  // e- check of the middle segment quality
113  // f- check if the distance between the reconstructed and the exstrapolated segments is ok
114 
115  // histo for efficiency with cuts a-/c-/d-
116  histos.push_back(ibooker.book2D(
117  "hEffGoodSegVsPosDen" + HistoName, "Eff vs local position (good) ", 25, -250., 250., 25, -250., 250.));
118  // histo for efficiency with cuts a-/b-/c-/d-/e-/f-
119  histos.push_back(ibooker.book2D("hEffGoodCloseSegVsPosNum" + HistoName,
120  "Eff vs local position (good and close segs) ",
121  25,
122  -250.,
123  250.,
124  25,
125  -250.,
126  250.));
127  if (detailedAnalysis) {
128  histos.push_back(
129  ibooker.book1D("hDistSegFromExtrap" + HistoName, "Distance segments from extrap position ", 200, 0., 200.));
130  // histo for efficiency from segment counting
131  histos.push_back(ibooker.book1D("hNaiveEffSeg" + HistoName, "Naive eff ", 10, 0., 10.));
132  // histo for efficiency with cuts a-/c-
133  histos.push_back(ibooker.book2D(
134  "hEffSegVsPosDen" + HistoName, "Eff vs local position (all) ", 25, -250., 250., 25, -250., 250.));
135  // histo for efficiency with cuts a-/b-/c-/d-
136  histos.push_back(
137  ibooker.book2D("hEffSegVsPosNum" + HistoName, "Eff vs local position ", 25, -250., 250., 25, -250., 250.));
138  // histo for efficiency with cuts a-/b-/c-/d-/e-
139  histos.push_back(ibooker.book2D(
140  "hEffGoodSegVsPosNum" + HistoName, "Eff vs local position (good segs) ", 25, -250., 250., 25, -250., 250.));
141  }
142  histosPerCh[chId] = histos;
143 }
144 
146  edm::LogVerbatim("DTDQM|DTMonitorModule|DTChamberEfficiencyTask")
147  << "[DTChamberEfficiencyTask] Analyze #Run: " << event.id().run() << " #Event: " << event.id().event();
148 
149  // Get the 4D rechit collection from the event
150  event.getByToken(recHits4DToken_, segs);
151 
152  int bottom = 0, top = 0;
153 
154  // Loop over all the chambers
155  vector<const DTChamber*>::const_iterator ch_it = dtGeom->chambers().begin();
156  vector<const DTChamber*>::const_iterator ch_end = dtGeom->chambers().end();
157  for (; ch_it != ch_end; ++ch_it) {
158  DTChamberId ch = (*ch_it)->id();
159  int wheel = ch.wheel();
160  int sector = ch.sector();
161  int station = ch.station();
162 
163  DTChamberId MidId(wheel, station, sector);
164 
165  // get efficiency for MB1 using MB2 and MB3
166  if (station == 1) {
167  bottom = 2;
168  top = 3;
169  }
170 
171  // get efficiency for MB2 using MB1 and MB3
172  if (station == 2) {
173  bottom = 1;
174  top = 3;
175  }
176 
177  // get efficiency for MB2 using MB2 and MB4
178  if (station == 3) {
179  bottom = 2;
180  top = 4;
181  }
182 
183  // get efficiency for MB4 using MB2 and MB3
184  if (station == 4) {
185  bottom = 2;
186  top = 3;
187  }
188 
189  // Select events with (good) segments in Bot and Top
190  DTChamberId BotId(wheel, bottom, sector);
191  DTChamberId TopId(wheel, top, sector);
192 
193  // Get segments in the bottom chambers (if any)
194  DTRecSegment4DCollection::range segsBot = segs->get(BotId);
195  int nSegsBot = segsBot.second - segsBot.first;
196  // check if any segments is there
197  if (nSegsBot == 0)
198  continue;
199 
200  vector<MonitorElement*> histos = histosPerCh[MidId];
201 
202  // Get segments in the top chambers (if any)
203  DTRecSegment4DCollection::range segsTop = segs->get(TopId);
204  int nSegsTop = segsTop.second - segsTop.first;
205 
206  // Select one segment for the bottom chamber
207  const DTRecSegment4D& bestBotSeg = getBestSegment(segsBot);
208 
209  // Select one segment for the top chamber
210  DTRecSegment4D* pBestTopSeg = nullptr;
211  if (nSegsTop > 0)
212  pBestTopSeg = const_cast<DTRecSegment4D*>(&getBestSegment(segsTop));
213  //if top chamber is MB4 sector 10, consider also sector 14
214  if (TopId.station() == 4 && TopId.sector() == 10) {
215  DTChamberId TopId14(wheel, top, 14);
216  DTRecSegment4DCollection::range segsTop14 = segs->get(TopId14);
217  int nSegsTop14 = segsTop14.second - segsTop14.first;
218  nSegsTop += nSegsTop14;
219  if (nSegsTop14) {
220  DTRecSegment4D* pBestTopSeg14 = const_cast<DTRecSegment4D*>(&getBestSegment(segsTop14));
221  // get best between sector 10 and 14
222  pBestTopSeg = const_cast<DTRecSegment4D*>(getBestSegment(pBestTopSeg, pBestTopSeg14));
223  }
224  }
225  if (!pBestTopSeg)
226  continue;
227  const DTRecSegment4D& bestTopSeg = *pBestTopSeg;
228 
229  DTRecSegment4DCollection::range segsMid = segs->get(MidId);
230  int nSegsMid = segsMid.second - segsMid.first;
231 
232  if (detailedAnalysis) {
233  // very trivial efficiency, just count segments
234  histos[3]->Fill(0);
235  if (nSegsMid > 0)
236  histos[3]->Fill(1);
237  }
238 
239  // get position at Mid by interpolating the position (not direction) of best
240  // segment in Bot and Top to Mid surface
241  LocalPoint posAtMid = interpolate(bestBotSeg, bestTopSeg, MidId);
242 
243  // is best segment good enough?
244  if (isGoodSegment(bestBotSeg) && isGoodSegment(bestTopSeg)) {
245  if (detailedAnalysis)
246  histos[4]->Fill(posAtMid.x(), posAtMid.y());
247  //check if interpolation fall inside middle chamber
248  if ((dtGeom->chamber(MidId))->surface().bounds().inside(posAtMid)) {
249  histos[0]->Fill(posAtMid.x(), posAtMid.y());
250  if (nSegsMid > 0) {
251  if (detailedAnalysis) {
252  histos[3]->Fill(2);
253  histos[5]->Fill(posAtMid.x(), posAtMid.y());
254  }
255 
256  const DTRecSegment4D& bestMidSeg = getBestSegment(segsMid);
257  // check if middle segments is good enough
258  if (isGoodSegment(bestMidSeg)) {
259  if (detailedAnalysis)
260  histos[6]->Fill(posAtMid.x(), posAtMid.y());
261  LocalPoint midSegPos = bestMidSeg.localPosition();
262 
263  // check if middle segments is also close enough
264  double dist;
265  if (bestMidSeg.hasPhi()) {
266  if (bestTopSeg.hasZed() && bestBotSeg.hasZed() && bestMidSeg.hasZed()) {
267  dist = (midSegPos - posAtMid).mag();
268  } else {
269  dist = fabs((midSegPos - posAtMid).x());
270  }
271  } else {
272  dist = fabs((midSegPos - posAtMid).y());
273  }
274  if (dist < theMinCloseDist) {
275  histos[1]->Fill(posAtMid.x(), posAtMid.y());
276  }
277  if (detailedAnalysis)
278  histos[2]->Fill(dist);
279  }
280  }
281  }
282  }
283  } // loop over stations
284 }
285 
286 // requirements : max number of hits and min chi2
289  unsigned int nHitBest = 0;
290  double chi2Best = 99999.;
291  for (DTRecSegment4DCollection::const_iterator seg = segs.first; seg != segs.second; ++seg) {
292  unsigned int nHits = ((*seg).hasPhi() ? (*seg).phiSegment()->recHits().size() : 0);
293  nHits += ((*seg).hasZed() ? (*seg).zSegment()->recHits().size() : 0);
294 
295  if (nHits == nHitBest) {
296  if ((*seg).chi2() / (*seg).degreesOfFreedom() < chi2Best) {
297  chi2Best = (*seg).chi2() / (*seg).degreesOfFreedom();
298  bestIter = seg;
299  }
300  } else if (nHits > nHitBest) {
301  nHitBest = nHits;
302  bestIter = seg;
303  }
304  }
305  return *bestIter;
306 }
307 
309  const DTRecSegment4D* s2) const {
310  if (!s1)
311  return s2;
312  if (!s2)
313  return s1;
314  unsigned int nHits1 = (s1->hasPhi() ? s1->phiSegment()->recHits().size() : 0);
315  nHits1 += (s1->hasZed() ? s1->zSegment()->recHits().size() : 0);
316 
317  unsigned int nHits2 = (s2->hasPhi() ? s2->phiSegment()->recHits().size() : 0);
318  nHits2 += (s2->hasZed() ? s2->zSegment()->recHits().size() : 0);
319 
320  if (nHits1 == nHits2) {
321  if (s1->chi2() / s1->degreesOfFreedom() < s2->chi2() / s2->degreesOfFreedom())
322  return s1;
323  else
324  return s2;
325  } else if (nHits1 > nHits2)
326  return s1;
327  return s2;
328 }
329 
331  const DTRecSegment4D& seg3,
332  const DTChamberId& id2) const {
333  // Get GlobalPoition of Seg in MB1
334  GlobalPoint gpos1 = (dtGeom->chamber(seg1.chamberId()))->toGlobal(seg1.localPosition());
335 
336  // Get GlobalPoition of Seg in MB3
337  GlobalPoint gpos3 = (dtGeom->chamber(seg3.chamberId()))->toGlobal(seg3.localPosition());
338 
339  // interpolate
340  // get all in MB2 frame
341  LocalPoint pos1 = (dtGeom->chamber(id2))->toLocal(gpos1);
342  LocalPoint pos3 = (dtGeom->chamber(id2))->toLocal(gpos3);
343 
344  // case 1: 1 and 3 has both projection. No problem
345 
346  // case 2: one projection is missing for one of the segments. Keep the other's segment position
347  if (!seg1.hasZed())
348  pos1 = LocalPoint(pos1.x(), pos3.y(), pos1.z());
349  if (!seg3.hasZed())
350  pos3 = LocalPoint(pos3.x(), pos1.y(), pos3.z());
351 
352  if (!seg1.hasPhi())
353  pos1 = LocalPoint(pos3.x(), pos1.y(), pos1.z());
354  if (!seg3.hasPhi())
355  pos3 = LocalPoint(pos1.x(), pos3.y(), pos3.z());
356 
357  // direction
358  LocalVector dir = (pos3 - pos1).unit(); // z points inward!
359  LocalPoint pos2 = pos1 + dir * pos1.z() / (-dir.z());
360 
361  return pos2;
362 }
363 
365  if (seg.chamberId().station() != 4 && !seg.hasZed())
366  return false;
367  unsigned int nHits = (seg.hasPhi() ? seg.phiSegment()->recHits().size() : 0);
368  nHits += (seg.hasZed() ? seg.zSegment()->recHits().size() : 0);
369  return (nHits >= theMinHitsSegment && seg.chi2() / seg.degreesOfFreedom() < theMinChi2NormSegment);
370 }
371 
372 // Local Variables:
373 // show-trailing-whitespace: t
374 // truncate-lines: t
375 // End:
LuminosityBlockID id() const
Log< level::Info, true > LogVerbatim
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< DTRecSegment4DCollection > recHits4DToken_
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:84
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:90
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
const DTChamberRecSegment2D * phiSegment() const
The superPhi segment: 0 if no phi projection available.
T y() const
Definition: PV3DBase.h:60
virtual DTChamberId chamberId() const
The (specific) DetId of the chamber on which the segment resides.
void dqmBeginRun(const edm::Run &run, const edm::EventSetup &setup) override
BeginRun.
LocalPoint localPosition() const override
Local position in Chamber frame.
edm::ESGetToken< DTGeometry, MuonGeometryRecord > muonGeomToken_
edm::Handle< DTRecSegment4DCollection > segs
bool getData(T &iHolder) const
Definition: EventSetup.h:128
~DTChamberEfficiencyTask() override
Destructor.
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
void beginLuminosityBlock(edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context) override
To reset the MEs.
T z() const
Definition: PV3DBase.h:61
void analyze(const edm::Event &event, const edm::EventSetup &setup) override
LocalPoint interpolate(const DTRecSegment4D &seg1, const DTRecSegment4D &seg3, const DTChamberId &MB2) const
Transition
Definition: Transition.h:12
std::vector< const TrackingRecHit * > recHits() const override
Access to component RecHits (if any)
bool hasPhi() const
Does it have the Phi projection?
DTChamberEfficiencyTask(const edm::ParameterSet &pset)
Constructor.
LocalVector toLocal(const reco::Track::Vector &v, const Surface &s)
bool hasZed() const
Does it have the Z projection?
const DTSLRecSegment2D * zSegment() const
The Z segment: 0 if not zed projection available.
bool isGoodSegment(const DTRecSegment4D &seg) const
std::map< DTChamberId, std::vector< MonitorElement * > > histosPerCh
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:177
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
LuminosityBlockNumber_t luminosityBlock() const
const DTRecSegment4D & getBestSegment(const DTRecSegment4DCollection::range &segs) const
caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple const cms::cuda::AtomicPairCounter GPUCACell const *__restrict__ uint32_t const *__restrict__ gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector const GPUCACell::OuterHitOfCell const int32_t nHits
int degreesOfFreedom() const override
Degrees of freedom of the segment fit.
std::string HistoName
int sector() const
Definition: DTChamberId.h:49
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
int station() const
Return the station number.
Definition: DTChamberId.h:42
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
void bookHistos(DQMStore::IBooker &ibooker, DTChamberId chId)
T x() const
Definition: PV3DBase.h:59
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
tuple size
Write out results.
Definition: Run.h:45
double chi2() const override
Chi2 of the segment fit.
Basic3DVector unit() const