CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
DTChamberEfficiencyTask Class Reference

#include <DTChamberEfficiencyTask.h>

Inheritance diagram for DTChamberEfficiencyTask:
one::DQMEDAnalyzer< edm::one::WatchLuminosityBlocks > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &setup) override
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context) override
 To reset the MEs. More...
 
void dqmBeginRun (const edm::Run &run, const edm::EventSetup &setup) override
 BeginRun. More...
 
 DTChamberEfficiencyTask (const edm::ParameterSet &pset)
 Constructor. More...
 
void endLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context) final
 
 ~DTChamberEfficiencyTask () override
 Destructor. More...
 
- Public Member Functions inherited from one::DQMEDAnalyzer< edm::one::WatchLuminosityBlocks >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Protected Member Functions

void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 

Private Member Functions

void bookHistos (DQMStore::IBooker &ibooker, DTChamberId chId)
 
const DTRecSegment4DgetBestSegment (const DTRecSegment4DCollection::range &segs) const
 
const DTRecSegment4DgetBestSegment (const DTRecSegment4D *s1, const DTRecSegment4D *s2) const
 
LocalPoint interpolate (const DTRecSegment4D &seg1, const DTRecSegment4D &seg3, const DTChamberId &MB2) const
 
bool isGoodSegment (const DTRecSegment4D &seg) const
 

Private Attributes

bool debug
 
bool detailedAnalysis
 
edm::ESHandle< DTGeometrydtGeom
 
std::map< DTChamberId, std::vector< MonitorElement * > > histosPerCh
 
bool onlineMonitor
 
edm::ParameterSet parameters
 
edm::EDGetTokenT< DTRecSegment4DCollectionrecHits4DToken_
 
edm::Handle< DTRecSegment4DCollectionsegs
 
double theMinChi2NormSegment
 
double theMinCloseDist
 
unsigned int theMinHitsSegment
 

Detailed Description

DQM Analysis of 4D DT segments, it produces plots about:

Class based on the code written by S. Lacaprara : RecoLocalMuon / DTSegment / test / DTEffAnalyzer.h

Author
G. Mila - INFN Torino

Definition at line 43 of file DTChamberEfficiencyTask.h.

Constructor & Destructor Documentation

DTChamberEfficiencyTask::DTChamberEfficiencyTask ( const edm::ParameterSet pset)

Constructor.

Definition at line 32 of file DTChamberEfficiencyTask.cc.

References debug, edm::ParameterSet::getUntrackedParameter(), and muonDTDigis_cfi::pset.

32  {
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 }
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< DTRecSegment4DCollection > recHits4DToken_
DTChamberEfficiencyTask::~DTChamberEfficiencyTask ( )
override

Destructor.

Definition at line 59 of file DTChamberEfficiencyTask.cc.

59  {
60 
61  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTChamberEfficiencyTask") << "[DTChamberEfficiencyTask] Destructor called!";
62 }

Member Function Documentation

void DTChamberEfficiencyTask::analyze ( const edm::Event event,
const edm::EventSetup setup 
)
override

Definition at line 153 of file DTChamberEfficiencyTask.cc.

References DTRecSegment4D::hasPhi(), DTRecSegment4D::hasZed(), plotFactory::histos, HiCaloJetParameters_cff::interpolate, DTRecSegment4D::localPosition(), mag(), DTChamberId::sector(), DTChamberId::station(), relativeConstraints::station, DTChamberId::wheel(), makeMuonMisalignmentScenario::wheel, PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by endLuminosityBlock().

153  {
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=nullptr;
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 }
edm::EDGetTokenT< DTRecSegment4DCollection > recHits4DToken_
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:102
LocalPoint localPosition() const override
Local position in Chamber frame.
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:117
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
T y() const
Definition: PV3DBase.h:63
edm::Handle< DTRecSegment4DCollection > segs
LocalPoint interpolate(const DTRecSegment4D &seg1, const DTRecSegment4D &seg3, const DTChamberId &MB2) const
bool hasPhi() const
Does it have the Phi projection?
bool hasZed() const
Does it have the Z projection?
edm::ESHandle< DTGeometry > dtGeom
bool isGoodSegment(const DTRecSegment4D &seg) const
std::map< DTChamberId, std::vector< MonitorElement * > > histosPerCh
const DTRecSegment4D & getBestSegment(const DTRecSegment4DCollection::range &segs) const
int sector() const
Definition: DTChamberId.h:61
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 DTChamberEfficiencyTask::beginLuminosityBlock ( edm::LuminosityBlock const &  lumiSeg,
edm::EventSetup const &  context 
)
override

To reset the MEs.

Definition at line 64 of file DTChamberEfficiencyTask.cc.

References trackerHits::histo, mps_fire::i, edm::LuminosityBlockBase::id(), edm::LuminosityBlockID::luminosityBlock(), genParticles_cff::map, and findQualityFiles::size.

64  {
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 }
size
Write out results.
std::map< DTChamberId, std::vector< MonitorElement * > > histosPerCh
void DTChamberEfficiencyTask::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &  context 
)
overrideprotected

Definition at line 89 of file DTChamberEfficiencyTask.cc.

References bookHistos(), and DQMStore::IBooker::setCurrentFolder().

Referenced by endLuminosityBlock().

89  {
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 }
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:102
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
edm::ESHandle< DTGeometry > dtGeom
void bookHistos(DQMStore::IBooker &ibooker, DTChamberId chId)
void DTChamberEfficiencyTask::bookHistos ( DQMStore::IBooker ibooker,
DTChamberId  chId 
)
private

Definition at line 104 of file DTChamberEfficiencyTask.cc.

References DQMStore::IBooker::book1D(), DQMStore::IBooker::book2D(), plotFactory::histos, DTChamberId::sector(), DQMStore::IBooker::setCurrentFolder(), DTChamberId::station(), relativeConstraints::station, DTChamberId::wheel(), and makeMuonMisalignmentScenario::wheel.

Referenced by endLuminosityBlock().

104  {
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 }
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
std::map< DTChamberId, std::vector< MonitorElement * > > histosPerCh
std::string HistoName
int sector() const
Definition: DTChamberId.h:61
int station() const
Return the station number.
Definition: DTChamberId.h:51
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
void DTChamberEfficiencyTask::dqmBeginRun ( const edm::Run run,
const edm::EventSetup setup 
)
override

BeginRun.

Definition at line 82 of file DTChamberEfficiencyTask.cc.

References edm::EventSetup::get().

82  {
83 
84  // Get the DT Geometry
85  setup.get<MuonGeometryRecord>().get(dtGeom);
86 
87 }
edm::ESHandle< DTGeometry > dtGeom
T get() const
Definition: EventSetup.h:71
void DTChamberEfficiencyTask::endLuminosityBlock ( edm::LuminosityBlock const &  lumiSeg,
edm::EventSetup const &  context 
)
inlinefinal
const DTRecSegment4D & DTChamberEfficiencyTask::getBestSegment ( const DTRecSegment4DCollection::range segs) const
private

Definition at line 302 of file DTChamberEfficiencyTask.cc.

References DTRecSegment4D::chi2().

Referenced by endLuminosityBlock().

302  {
304  unsigned int nHitBest=0;
305  double chi2Best=99999.;
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 }
edm::Handle< DTRecSegment4DCollection > segs
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
const DTRecSegment4D * DTChamberEfficiencyTask::getBestSegment ( const DTRecSegment4D s1,
const DTRecSegment4D s2 
) const
private

Definition at line 325 of file DTChamberEfficiencyTask.cc.

References DTRecSegment4D::chi2(), DTRecSegment4D::degreesOfFreedom(), DTRecSegment4D::hasPhi(), DTRecSegment4D::hasZed(), DTRecSegment4D::phiSegment(), DTRecSegment2D::recHits(), indexGen::s2, and DTRecSegment4D::zSegment().

326  {
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 }
const DTChamberRecSegment2D * phiSegment() const
The superPhi segment: 0 if no phi projection available.
double chi2() const override
Chi2 of the segment fit.
bool hasPhi() const
Does it have the Phi projection?
bool hasZed() const
Does it have the Z projection?
const DTSLRecSegment2D * zSegment() const
The Z segment: 0 if not zed projection available.
std::vector< const TrackingRecHit * > recHits() const override
Access to component RecHits (if any)
int degreesOfFreedom() const override
Degrees of freedom of the segment fit.
LocalPoint DTChamberEfficiencyTask::interpolate ( const DTRecSegment4D seg1,
const DTRecSegment4D seg3,
const DTChamberId MB2 
) const
private

Definition at line 347 of file DTChamberEfficiencyTask.cc.

References DTRecSegment4D::chamberId(), dir, DTRecSegment4D::hasPhi(), DTRecSegment4D::hasZed(), DTRecSegment4D::localPosition(), toLocal(), csvLumiCalc::unit, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by endLuminosityBlock().

349  {
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 }
LocalPoint localPosition() const override
Local position in Chamber frame.
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:32
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:117
T y() const
Definition: PV3DBase.h:63
virtual DTChamberId chamberId() const
The (specific) DetId of the chamber on which the segment resides.
LocalVector toLocal(const reco::Track::Vector &v, const Surface &s)
T z() const
Definition: PV3DBase.h:64
bool hasPhi() const
Does it have the Phi projection?
bool hasZed() const
Does it have the Z projection?
edm::ESHandle< DTGeometry > dtGeom
dbl *** dir
Definition: mlp_gen.cc:35
T x() const
Definition: PV3DBase.h:62
bool DTChamberEfficiencyTask::isGoodSegment ( const DTRecSegment4D seg) const
private

Definition at line 378 of file DTChamberEfficiencyTask.cc.

References DTRecSegment4D::chamberId(), DTRecSegment4D::chi2(), DTRecSegment4D::degreesOfFreedom(), DTRecSegment4D::hasPhi(), DTRecSegment4D::hasZed(), DTRecSegment4D::phiSegment(), DTRecSegment2D::recHits(), DTChamberId::station(), and DTRecSegment4D::zSegment().

Referenced by endLuminosityBlock().

378  {
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 &&
384 }
const DTChamberRecSegment2D * phiSegment() const
The superPhi segment: 0 if no phi projection available.
double chi2() const override
Chi2 of the segment fit.
virtual DTChamberId chamberId() const
The (specific) DetId of the chamber on which the segment resides.
bool hasPhi() const
Does it have the Phi projection?
bool hasZed() const
Does it have the Z projection?
const DTSLRecSegment2D * zSegment() const
The Z segment: 0 if not zed projection available.
std::vector< const TrackingRecHit * > recHits() const override
Access to component RecHits (if any)
int station() const
Return the station number.
Definition: DTChamberId.h:51
int degreesOfFreedom() const override
Degrees of freedom of the segment fit.

Member Data Documentation

bool DTChamberEfficiencyTask::debug
private
bool DTChamberEfficiencyTask::detailedAnalysis
private

Definition at line 82 of file DTChamberEfficiencyTask.h.

edm::ESHandle<DTGeometry> DTChamberEfficiencyTask::dtGeom
private

Definition at line 95 of file DTChamberEfficiencyTask.h.

std::map<DTChamberId, std::vector<MonitorElement*> > DTChamberEfficiencyTask::histosPerCh
private

Definition at line 89 of file DTChamberEfficiencyTask.h.

bool DTChamberEfficiencyTask::onlineMonitor
private

Definition at line 80 of file DTChamberEfficiencyTask.h.

edm::ParameterSet DTChamberEfficiencyTask::parameters
private
edm::EDGetTokenT<DTRecSegment4DCollection> DTChamberEfficiencyTask::recHits4DToken_
private

Definition at line 85 of file DTChamberEfficiencyTask.h.

edm::Handle<DTRecSegment4DCollection> DTChamberEfficiencyTask::segs
private

Definition at line 96 of file DTChamberEfficiencyTask.h.

Referenced by endLuminosityBlock().

double DTChamberEfficiencyTask::theMinChi2NormSegment
private

Definition at line 92 of file DTChamberEfficiencyTask.h.

double DTChamberEfficiencyTask::theMinCloseDist
private

Definition at line 93 of file DTChamberEfficiencyTask.h.

unsigned int DTChamberEfficiencyTask::theMinHitsSegment
private

Definition at line 91 of file DTChamberEfficiencyTask.h.