CMS 3D CMS Logo

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