CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTChamberEfficiencyTask.cc
Go to the documentation of this file.
1 
2 
3 /*
4  * See header file for a description of this class.
5  *
6  * \author G. Mila - INFN Torino
7  */
8 
9 
10 
12 
13 //Framework
17 
21 
22 //Geometry
24 
25 
26 #include <iterator>
27 #include <iostream>
28 #include <cmath>
29 using namespace edm;
30 using namespace std;
31 
32 
33 
35 
36  debug = pset.getUntrackedParameter<bool>("debug",false);
37 
38  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTChamberEfficiencyTask") << "[DTChamberEfficiencyTask] Constructor called!";
39 
40  // Get the DQM needed services
41  theDbe = edm::Service<DQMStore>().operator->();
42 
43  theDbe->setCurrentFolder("DT/DTChamberEfficiencyTask");
44 
45  parameters = pset;
46 
47 }
48 
49 
51 
52  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTChamberEfficiencyTask") << "[DTChamberEfficiencyTask] Destructor called!";
53 }
54 
55 
57 
58  // the name of the 4D rec hits collection
59  recHits4DToken_ = consumes<DTRecSegment4DCollection>(
60  edm::InputTag(parameters.getParameter<string>("recHits4DLabel")));
61 
62  // parameters to use for the segment quality check
63  theMinHitsSegment = static_cast<unsigned int>(parameters.getParameter<int>("minHitsSegment"));
64  theMinChi2NormSegment = parameters.getParameter<double>("minChi2NormSegment");
65  // parameter to use for the exstrapolated segment check
66  theMinCloseDist = parameters.getParameter<double>("minCloseDist");
67 
68  // the running modality
69  onlineMonitor = parameters.getUntrackedParameter<bool>("onlineMonitor");
70 
71  // the analysis mode
72  detailedAnalysis = parameters.getUntrackedParameter<bool>("detailedAnalysis");
73 
74 }
75 
76 
78 
79  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTChamberEfficiencyTask")<<"[DTChamberEfficiencyTask]: Begin of LS transition";
80 
81  if(lumiSeg.id().luminosityBlock()%parameters.getUntrackedParameter<int>("ResetCycle", 3) == 0 && onlineMonitor) {
82  for(map<DTChamberId, vector<MonitorElement*> > ::const_iterator histo = histosPerCh.begin();
83  histo != histosPerCh.end();
84  histo++) {
85  int size = (*histo).second.size();
86  for(int i=0; i<size; i++){
87  (*histo).second[i]->Reset();
88  }
89  }
90  }
91 
92 }
93 
94 
96 
97  // Get the DT Geometry
98  setup.get<MuonGeometryRecord>().get(dtGeom);
99 
100  // Loop over all the chambers
101  vector<const DTChamber*>::const_iterator ch_it = dtGeom->chambers().begin();
102  vector<const DTChamber*>::const_iterator ch_end = dtGeom->chambers().end();
103  for (; ch_it != ch_end; ++ch_it) {
104  // histo booking
105  bookHistos((*ch_it)->id());
106  }
107 
108 }
109 
110 
111 
113 
114  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTChamberEfficiencyTask")<<"[DTChamberEfficiencyTask] endjob called!";
115 }
116 
117 
118 
119 // Book a set of histograms for a given Layer
121 
122  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTChamberEfficiencyTask") << " Booking histos for CH : " << chId;
123 
124  // Compose the chamber name
125  stringstream wheel; wheel << chId.wheel();
126  stringstream station; station << chId.station();
127  stringstream sector; sector << chId.sector();
128 
129  string HistoName =
130  "_W" + wheel.str() +
131  "_St" + station.str() +
132  "_Sec" + sector.str();
133 
134  theDbe->setCurrentFolder("DT/01-DTChamberEfficiency/Task/Wheel" + wheel.str() +
135  "/Sector" + sector.str() +
136  "/Station" + station.str());
137 
138  // Create the monitor elements
139  vector<MonitorElement *> histos;
140 
141  //efficiency selection cuts
142  // a- number of segments of the top chamber > 0 && number of segments of the bottom chamber > 0
143  // b- number of segments of the middle chamber > 0
144  // c- check of the top and bottom segment quality
145  // d- check if interpolation falls inside the middle chamber
146  // e- check of the middle segment quality
147  // f- check if the distance between the reconstructed and the exstrapolated segments is ok
148 
149 
150  // histo for efficiency with cuts a-/c-/d-
151  histos.push_back(theDbe->book2D("hEffGoodSegVsPosDen"+HistoName,"Eff vs local position (good) ",25,-250.,250., 25,-250.,250.));
152  // histo for efficiency with cuts a-/b-/c-/d-/e-/f-
153  histos.push_back(theDbe->book2D("hEffGoodCloseSegVsPosNum"+HistoName, "Eff vs local position (good and close segs) ", 25,-250.,250., 25,-250.,250.));
154  if(detailedAnalysis){
155  histos.push_back(theDbe->book1D("hDistSegFromExtrap"+HistoName, "Distance segments from extrap position ",200,0.,200.));
156  // histo for efficiency from segment counting
157  histos.push_back(theDbe->book1D("hNaiveEffSeg"+HistoName, "Naive eff ",10,0.,10.));
158  // histo for efficiency with cuts a-/c-
159  histos.push_back(theDbe->book2D("hEffSegVsPosDen"+HistoName,"Eff vs local position (all) ",25,-250.,250., 25,-250.,250.));
160  // histo for efficiency with cuts a-/b-/c-/d-
161  histos.push_back(theDbe->book2D("hEffSegVsPosNum"+HistoName, "Eff vs local position ",25,-250.,250., 25,-250.,250.));
162  // histo for efficiency with cuts a-/b-/c-/d-/e-
163  histos.push_back(theDbe->book2D("hEffGoodSegVsPosNum"+HistoName, "Eff vs local position (good segs) ", 25,-250.,250., 25,-250.,250.));
164  }
165  histosPerCh[chId] = histos;
166 }
167 
168 
170 
171  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTChamberEfficiencyTask") << "[DTChamberEfficiencyTask] Analyze #Run: " << event.id().run()
172  << " #Event: " << event.id().event();
173 
174  // Get the 4D rechit collection from the event
175  event.getByToken(recHits4DToken_, segs);
176 
177  int bottom=0, top=0;
178 
179 
180  // Loop over all the chambers
181  vector<const DTChamber*>::const_iterator ch_it = dtGeom->chambers().begin();
182  vector<const DTChamber*>::const_iterator ch_end = dtGeom->chambers().end();
183  for (; ch_it != ch_end; ++ch_it) {
184 
185  DTChamberId ch = (*ch_it)->id();
186  int wheel = ch.wheel();
187  int sector = ch.sector();
188  int station = ch.station();
189 
190 
191  DTChamberId MidId(wheel, station, sector);
192 
193  // get efficiency for MB1 using MB2 and MB3
194  if( station == 1 ) {
195  bottom = 2;
196  top = 3;
197  }
198 
199  // get efficiency for MB2 using MB1 and MB3
200  if( station == 2 ) {
201  bottom = 1;
202  top = 3;
203  }
204 
205  // get efficiency for MB2 using MB2 and MB4
206  if( station == 3 ) {
207  bottom = 2;
208  top = 4;
209  }
210 
211  // get efficiency for MB4 using MB2 and MB3
212  if( station == 4 ) {
213  bottom = 2;
214  top = 3;
215  }
216 
217  // Select events with (good) segments in Bot and Top
218  DTChamberId BotId(wheel, bottom, sector);
219  DTChamberId TopId(wheel, top, sector);
220 
221  // Get segments in the bottom chambers (if any)
222  DTRecSegment4DCollection::range segsBot= segs->get(BotId);
223  int nSegsBot=segsBot.second-segsBot.first;
224  // check if any segments is there
225  if (nSegsBot==0) continue;
226 
227  vector<MonitorElement *> histos = histosPerCh[MidId];
228 
229  // Get segments in the top chambers (if any)
230  DTRecSegment4DCollection::range segsTop= segs->get(TopId);
231  int nSegsTop=segsTop.second-segsTop.first;
232 
233  // Select one segment for the bottom chamber
234  const DTRecSegment4D& bestBotSeg= getBestSegment(segsBot);
235 
236  // Select one segment for the top chamber
237  DTRecSegment4D* pBestTopSeg=0;
238  if (nSegsTop>0)
239  pBestTopSeg = const_cast<DTRecSegment4D*>(&getBestSegment(segsTop));
240  //if top chamber is MB4 sector 10, consider also sector 14
241  if (TopId.station() == 4 && TopId.sector() == 10) {
242  DTChamberId TopId14(wheel, top, 14);
243  DTRecSegment4DCollection::range segsTop14= segs->get(TopId14);
244  int nSegsTop14=segsTop14.second-segsTop14.first;
245  nSegsTop+=nSegsTop14;
246  if (nSegsTop14) {
247  DTRecSegment4D* pBestTopSeg14 = const_cast<DTRecSegment4D*>(&getBestSegment(segsTop14));
248  // get best between sector 10 and 14
249  pBestTopSeg = const_cast<DTRecSegment4D*>(getBestSegment(pBestTopSeg, pBestTopSeg14));
250  }
251  }
252  if (!pBestTopSeg) continue;
253  const DTRecSegment4D& bestTopSeg= *pBestTopSeg;
254 
255  DTRecSegment4DCollection::range segsMid= segs->get(MidId);
256  int nSegsMid=segsMid.second-segsMid.first;
257 
258  if(detailedAnalysis){
259  // very trivial efficiency, just count segments
260  histos[3]->Fill(0);
261  if (nSegsMid>0) histos[3]->Fill(1);
262  }
263 
264  // get position at Mid by interpolating the position (not direction) of best
265  // segment in Bot and Top to Mid surface
266  LocalPoint posAtMid = interpolate(bestBotSeg, bestTopSeg, MidId);
267 
268  // is best segment good enough?
269  if (isGoodSegment(bestBotSeg) && isGoodSegment(bestTopSeg)) {
270  if(detailedAnalysis)
271  histos[4]->Fill(posAtMid.x(),posAtMid.y());
272  //check if interpolation fall inside middle chamber
273  if ((dtGeom->chamber(MidId))->surface().bounds().inside(posAtMid)) {
274  histos[0]->Fill(posAtMid.x(),posAtMid.y());
275  if (nSegsMid>0) {
276 
277  if(detailedAnalysis){
278  histos[3]->Fill(2);
279  histos[5]->Fill(posAtMid.x(),posAtMid.y());
280  }
281 
282  const DTRecSegment4D& bestMidSeg= getBestSegment(segsMid);
283  // check if middle segments is good enough
284  if (isGoodSegment(bestMidSeg)) {
285 
286  if(detailedAnalysis)
287  histos[6]->Fill(posAtMid.x(),posAtMid.y());
288  LocalPoint midSegPos=bestMidSeg.localPosition();
289 
290  // check if middle segments is also close enough
291  double dist;
292  if (bestMidSeg.hasPhi()) {
293  if (bestTopSeg.hasZed() && bestBotSeg.hasZed() && bestMidSeg.hasZed()) {
294  dist = (midSegPos-posAtMid).mag();
295  } else {
296  dist = fabs((midSegPos-posAtMid).x());
297  }
298  } else {
299  dist = fabs((midSegPos-posAtMid).y());
300  }
301  if (dist < theMinCloseDist ) {
302  histos[1]->Fill(posAtMid.x(),posAtMid.y());
303  }
304  if(detailedAnalysis)
305  histos[2]->Fill(dist);
306  }
307  }
308  }
309  }
310  }// loop over stations
311 
312 }
313 
314 
315 
316 
317 // requirements : max number of hits and min chi2
320  unsigned int nHitBest=0;
321  double chi2Best=99999.;
322  for (DTRecSegment4DCollection::const_iterator seg=segs.first ;
323  seg!=segs.second ; ++seg ) {
324  unsigned int nHits= ((*seg).hasPhi() ? (*seg).phiSegment()->recHits().size() : 0 ) ;
325  nHits+= ((*seg).hasZed() ? (*seg).zSegment()->recHits().size() : 0 );
326 
327  if (nHits==nHitBest) {
328  if ((*seg).chi2()/(*seg).degreesOfFreedom() < chi2Best ) {
329  chi2Best=(*seg).chi2()/(*seg).degreesOfFreedom();
330  bestIter = seg;
331  }
332  }
333  else if (nHits>nHitBest) {
334  nHitBest=nHits;
335  bestIter = seg;
336  }
337  }
338  return *bestIter;
339 }
340 
342  const DTRecSegment4D* s2) const{
343 
344  if (!s1) return s2;
345  if (!s2) return s1;
346  unsigned int nHits1= (s1->hasPhi() ? s1->phiSegment()->recHits().size() : 0 ) ;
347  nHits1+= (s1->hasZed() ? s1->zSegment()->recHits().size() : 0 );
348 
349  unsigned int nHits2= (s2->hasPhi() ? s2->phiSegment()->recHits().size() : 0 ) ;
350  nHits2+= (s2->hasZed() ? s2->zSegment()->recHits().size() : 0 );
351 
352  if (nHits1==nHits2) {
353  if (s1->chi2()/s1->degreesOfFreedom() < s2->chi2()/s2->degreesOfFreedom() )
354  return s1;
355  else
356  return s2;
357  }
358  else if (nHits1>nHits2) return s1;
359  return s2;
360 }
361 
362 
364  const DTRecSegment4D& seg3,
365  const DTChamberId& id2) const {
366  // Get GlobalPoition of Seg in MB1
367  GlobalPoint gpos1=(dtGeom->chamber(seg1.chamberId()))->toGlobal(seg1.localPosition());
368 
369  // Get GlobalPoition of Seg in MB3
370  GlobalPoint gpos3=(dtGeom->chamber(seg3.chamberId()))->toGlobal(seg3.localPosition());
371 
372  // interpolate
373  // get all in MB2 frame
374  LocalPoint pos1=(dtGeom->chamber(id2))->toLocal(gpos1);
375  LocalPoint pos3=(dtGeom->chamber(id2))->toLocal(gpos3);
376 
377  // case 1: 1 and 3 has both projection. No problem
378 
379  // case 2: one projection is missing for one of the segments. Keep the other's segment position
380  if (!seg1.hasZed()) pos1=LocalPoint(pos1.x(),pos3.y(),pos1.z());
381  if (!seg3.hasZed()) pos3=LocalPoint(pos3.x(),pos1.y(),pos3.z());
382 
383  if (!seg1.hasPhi()) pos1=LocalPoint(pos3.x(),pos1.y(),pos1.z());
384  if (!seg3.hasPhi()) pos3=LocalPoint(pos1.x(),pos3.y(),pos3.z());
385 
386  // direction
387  LocalVector dir = (pos3-pos1).unit(); // z points inward!
388  LocalPoint pos2 = pos1+dir*pos1.z()/(-dir.z());
389 
390  return pos2;
391 }
392 
393 
395  if (seg.chamberId().station()!=4 && !seg.hasZed()) return false;
396  unsigned int nHits= (seg.hasPhi() ? seg.phiSegment()->recHits().size() : 0 ) ;
397  nHits+= (seg.hasZed() ? seg.zSegment()->recHits().size() : 0 );
398  return ( nHits >= theMinHitsSegment &&
399  seg.chi2()/seg.degreesOfFreedom() < theMinChi2NormSegment );
400 }
401 
402 // Local Variables:
403 // show-trailing-whitespace: t
404 // truncate-lines: t
405 // End:
LuminosityBlockID id() const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
dictionary parameters
Definition: Parameters.py:2
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
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
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.
LocalVector toLocal(const reco::Track::Vector &v, const Surface &s)
tuple s2
Definition: indexGen.py:106
void bookHistos()
Definition: Histogram.h:33
void bookHistos(DTChamberId chId)
void analyze(const edm::Event &event, const edm::EventSetup &setup)
virtual std::vector< const TrackingRecHit * > recHits() const
Access to component RecHits (if any)
string unit
Definition: csvLumiCalc.py:46
T z() const
Definition: PV3DBase.h:64
LocalPoint interpolate(const DTRecSegment4D &seg1, const DTRecSegment4D &seg3, const DTChamberId &MB2) const
bool hasPhi() const
Does it have the Phi projection?
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
virtual LocalPoint localPosition() const
Local position in Chamber frame.
DTChamberEfficiencyTask(const edm::ParameterSet &pset)
Constructor.
bool hasZed() const
Does it have the Z projection?
string top
Definition: fff_deleter.py:272
#define debug
Definition: HDRShower.cc:19
const DTSLRecSegment2D * zSegment() const
The Z segment: 0 if not zed projection available.
bool isGoodSegment(const DTRecSegment4D &seg) const
const T & get() const
Definition: EventSetup.h:55
LuminosityBlockNumber_t luminosityBlock() const
const DTRecSegment4D & getBestSegment(const DTRecSegment4DCollection::range &segs) const
std::string HistoName
int sector() const
Definition: DTChamberId.h:61
void beginRun(const edm::Run &run, const edm::EventSetup &setup)
BeginRun.
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
virtual ~DTChamberEfficiencyTask()
Destructor.
dbl *** dir
Definition: mlp_gen.cc:35
Definition: DDAxes.h:10
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
T x() const
Definition: PV3DBase.h:62
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
tuple size
Write out results.
Definition: Run.h:41
void beginLuminosityBlock(edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context)
To reset the MEs.