CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
DTSegment4DQuality Class Reference

#include <DTSegment4DQuality.h>

Inheritance diagram for DTSegment4DQuality:
edm::EDAnalyzer

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup)
 Perform the real analysis. More...
 
 DTSegment4DQuality (const edm::ParameterSet &pset)
 Constructor. More...
 
void endJob ()
 
void endLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &c)
 
virtual ~DTSegment4DQuality ()
 Destructor. More...
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Attributes

DQMStoredbe_
 
bool debug
 
bool doall
 
HRes4DHith4DHit
 
HRes4DHith4DHit_W0
 
HRes4DHith4DHit_W1
 
HRes4DHith4DHit_W2
 
HRes4DHith4DHitWS [3][4]
 
HEff4DHithEff_All
 
HEff4DHithEff_W0
 
HEff4DHithEff_W1
 
HEff4DHithEff_W2
 
HEff4DHithEffWS [3][4]
 
MonitorElementhHitMult [3][4]
 
MonitorElementht0 [3][4]
 
bool local
 
std::string rootFileName
 
edm::InputTag segment4DLabel
 
double sigmaResAlpha
 
double sigmaResBeta
 
double sigmaResX
 
double sigmaResY
 
edm::InputTag simHitLabel
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Basic analyzer class which accesses 4D DTSegments and plots resolution comparing reconstructed and simulated quantities

Only true 4D segments are considered. Station 4 segments are not looked at. FIXME: Add flag to consider also segments with only phi view? Possible bias?

Residual/pull plots are filled for the reco segment with alpha closest to the simulated muon direction (defined from muon simhits in the chamber).

Efficiencies are defined as reconstructed 4D segments with alpha, beta, x, y, within 5 sigma relative to the sim muon, with sigmas specified in the config. Note that loss of even only one of the two views is considered as inefficiency!

Author
S. Bolognesi and G. Cerminara - INFN Torino

Definition at line 43 of file DTSegment4DQuality.h.

Constructor & Destructor Documentation

DTSegment4DQuality::DTSegment4DQuality ( const edm::ParameterSet pset)

Constructor.

Definition at line 39 of file DTSegment4DQuality.cc.

References dbe_, debug, DTHitQualityUtils::debug, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), mergeVDriftHistosByStation::name, cppFunctionSkipper::operator, dtTPAnalyzer_cfg::rootFileName, alignCSCRings::s, and w().

39  {
40  // Get the debug parameter for verbose output
41  debug = pset.getUntrackedParameter<bool>("debug");
43 
44  rootFileName = pset.getUntrackedParameter<string>("rootFileName");
45  // the name of the simhit collection
46  simHitLabel = pset.getUntrackedParameter<InputTag>("simHitLabel");
47  // the name of the 4D rec hit collection
48  segment4DLabel = pset.getUntrackedParameter<InputTag>("segment4DLabel");
49 
50  //sigma resolution on position
51  sigmaResX = pset.getParameter<double>("sigmaResX");
52  sigmaResY = pset.getParameter<double>("sigmaResY");
53  //sigma resolution on angle
54  sigmaResAlpha = pset.getParameter<double>("sigmaResAlpha");
55  sigmaResBeta = pset.getParameter<double>("sigmaResBeta");
56  doall = pset.getUntrackedParameter<bool>("doall", false);
57  local = pset.getUntrackedParameter<bool>("local", false);
58 
59 
60  // Create the root file
61  //theFile = new TFile(rootFileName.c_str(), "RECREATE");
62  //theFile->cd();
63 // ----------------------
64  // get hold of back-end interface
65  dbe_ = 0;
66  dbe_ = Service<DQMStore>().operator->();
67  if ( dbe_ ) {
68  if (debug) {
69  dbe_->setVerbose(1);
70  } else {
71  dbe_->setVerbose(0);
72  }
73  }
74  if ( dbe_ ) {
75  if ( debug ) dbe_->showDirStructure();
76  }
77 
78  h4DHit= new HRes4DHit ("All",dbe_,doall,local);
79  h4DHit_W0= new HRes4DHit ("W0",dbe_,doall,local);
80  h4DHit_W1= new HRes4DHit ("W1",dbe_,doall,local);
81  h4DHit_W2= new HRes4DHit ("W2",dbe_,doall,local);
82 
83  if(doall) {
84  hEff_All= new HEff4DHit ("All",dbe_);
85  hEff_W0= new HEff4DHit ("W0",dbe_);
86  hEff_W1= new HEff4DHit ("W1",dbe_);
87  hEff_W2= new HEff4DHit ("W2",dbe_);
88  }
89 
90  if (local) {
91  // Plots with finer granularity, not to be included in DQM
92  TString name="W";
93  for (long w=0;w<=2;++w) {
94  for (long s=1;s<=4;++s){
95  // FIXME station 4 is not filled
96  TString nameWS=(name+w+"_St"+s);
97  h4DHitWS[w][s-1] = new HRes4DHit(nameWS.Data(),dbe_,doall,local);
98  hEffWS[w][s-1] = new HEff4DHit (nameWS.Data(),dbe_);
99  dbe_->setCurrentFolder("DT/4DSegments/");
100  hHitMult[w][s-1] = dbe_->book2D("4D_" +nameWS+ "_hNHits", "NHits",12,0,12, 6,0,6);
101  ht0[w][s-1] = dbe_->book2D("4D_" +nameWS+ "_ht0", "t0",200,-25,25,200,-25,25);
102  }
103  }
104  }
105 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::InputTag segment4DLabel
HRes4DHit * h4DHitWS[3][4]
MonitorElement * hHitMult[3][4]
A set of histograms for efficiency 4D RecHits.
Definition: Histograms.h:1116
void setVerbose(unsigned level)
Definition: DQMStore.cc:393
edm::InputTag simHitLabel
HEff4DHit * hEffWS[3][4]
void showDirStructure(void) const
Definition: DQMStore.cc:2761
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:845
MonitorElement * ht0[3][4]
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
T w() const
DTSegment4DQuality::~DTSegment4DQuality ( )
virtual

Destructor.

Definition at line 108 of file DTSegment4DQuality.cc.

108  {
109 
110 }

Member Function Documentation

void DTSegment4DQuality::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
)
virtual

Perform the real analysis.

Implements edm::EDAnalyzer.

Definition at line 142 of file DTSegment4DQuality.cc.

References abs, funct::cos(), gather_cfg::cout, debug, PV3DBase< T, PVType, FrameType >::eta(), HRes4DHit::Fill(), HEff4DHit::Fill(), DTHitQualityUtils::findMuSimSegment(), DTHitQualityUtils::findMuSimSegmentDirAndPos(), DTHitQualityUtils::findSegmentAlphaAndBeta(), edm::EventSetup::get(), interpolateCardsSimple::histo, edm::HandleBase::isValid(), DTRecSegment4D::localDirection(), DTRecSegment2D::localDirection(), DTRecSegment4D::localDirectionError(), DTRecSegment2D::localDirectionError(), DTRecSegment4D::localPosition(), DTRecSegment2D::localPosition(), DTRecSegment4D::localPositionError(), DTRecSegment2D::localPositionError(), DTHitQualityUtils::mapMuSimHitsPerWire(), DTHitQualityUtils::mapSimHitsPerWire(), PV3DBase< T, PVType, FrameType >::phi(), DTRecSegment2D::recHits(), trackerHits::simHits, mathSSE::sqrt(), DTSLRecSegment2D::superLayerId(), DTRecSegment2D::t0(), PV3DBase< T, PVType, FrameType >::theta(), GeomDet::toGlobal(), GeomDet::toLocal(), PV3DBase< T, PVType, FrameType >::x(), LocalError::xx(), PV3DBase< T, PVType, FrameType >::y(), LocalError::yy(), PV3DBase< T, PVType, FrameType >::z(), and DTRecSegment4D::zSegment().

142  {
143  //theFile->cd();
144 
145  // Get the DT Geometry
146  ESHandle<DTGeometry> dtGeom;
147  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
148 
149  // Get the SimHit collection from the event
151  event.getByLabel(simHitLabel, simHits); //FIXME: second string to be removed
152 
153  //Map simHits by chamber
154  map<DTChamberId, PSimHitContainer > simHitsPerCh;
155  for(PSimHitContainer::const_iterator simHit = simHits->begin();
156  simHit != simHits->end(); simHit++){
157  // Create the id of the chamber (the simHits in the DT known their wireId)
158  DTChamberId chamberId = (((DTWireId(simHit->detUnitId())).layerId()).superlayerId()).chamberId();
159  // Fill the map
160  simHitsPerCh[chamberId].push_back(*simHit);
161  }
162 
163  // Get the 4D rechits from the event
165  event.getByLabel(segment4DLabel, segment4Ds);
166 
167  if(!segment4Ds.isValid()) {
168  if(debug) cout << "[DTSegment4DQuality]**Warning: no 4D Segments with label: " <<segment4DLabel
169  << " in this event, skipping!" << endl;
170  return;
171  }
172 
173  // Loop over all chambers containing a segment
175  for (chamberId = segment4Ds->id_begin();
176  chamberId != segment4Ds->id_end();
177  ++chamberId){
178 
179  if((*chamberId).station() == 4)
180  continue; //use DTSegment2DSLPhiQuality to analyze MB4 performaces
181 
182  //------------------------- simHits ---------------------------//
183  //Get simHits of each chamber
184  PSimHitContainer simHits = simHitsPerCh[(*chamberId)];
185 
186  // Map simhits per wire
187  map<DTWireId, PSimHitContainer > simHitsPerWire = DTHitQualityUtils::mapSimHitsPerWire(simHits);
188  map<DTWireId, const PSimHit*> muSimHitPerWire = DTHitQualityUtils::mapMuSimHitsPerWire(simHitsPerWire);
189  int nMuSimHit = muSimHitPerWire.size();
190  if(nMuSimHit == 0 || nMuSimHit == 1) {
191  if(debug && nMuSimHit == 1)
192  cout << "[DTSegment4DQuality] Only " << nMuSimHit << " mu SimHit in this chamber, skipping!" << endl;
193  continue; // If no or only one mu SimHit is found skip this chamber
194  }
195  if(debug)
196  cout << "=== Chamber " << (*chamberId) << " has " << nMuSimHit << " SimHits" << endl;
197 
198  //Find outer and inner mu SimHit to build a segment
199  pair<const PSimHit*, const PSimHit*> inAndOutSimHit = DTHitQualityUtils::findMuSimSegment(muSimHitPerWire);
200 
201  //Find direction and position of the sim Segment in Chamber RF
202  pair<LocalVector, LocalPoint> dirAndPosSimSegm = DTHitQualityUtils::findMuSimSegmentDirAndPos(inAndOutSimHit,
203  (*chamberId),&(*dtGeom));
204 
205  LocalVector simSegmLocalDir = dirAndPosSimSegm.first;
206  LocalPoint simSegmLocalPos = dirAndPosSimSegm.second;
207  const DTChamber* chamber = dtGeom->chamber(*chamberId);
208  GlobalPoint simSegmGlobalPos = chamber->toGlobal(simSegmLocalPos);
209  GlobalVector simSegmGlobalDir = chamber->toGlobal(simSegmLocalDir);
210 
211  //phi and theta angle of simulated segment in Chamber RF
212  float alphaSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).first;
213  float betaSimSeg = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegmLocalDir).second;
214  //x,y position of simulated segment in Chamber RF
215  float xSimSeg = simSegmLocalPos.x();
216  float ySimSeg = simSegmLocalPos.y();
217  //Position (in eta,phi coordinates) in lobal RF
218  float etaSimSeg = simSegmGlobalPos.eta();
219  float phiSimSeg = simSegmGlobalPos.phi();
220 
221  if(debug)
222  cout<<" Simulated segment: local direction "<<simSegmLocalDir<<endl
223  <<" local position "<<simSegmLocalPos<<endl
224  <<" alpha "<<alphaSimSeg<<endl
225  <<" beta "<<betaSimSeg<<endl;
226 
227  //---------------------------- recHits --------------------------//
228  // Get the range of rechit for the corresponding chamberId
229  bool recHitFound = false;
230  DTRecSegment4DCollection::range range = segment4Ds->get(*chamberId);
231  int nsegm = distance(range.first, range.second);
232  if(debug)
233  cout << " Chamber: " << *chamberId << " has " << nsegm
234  << " 4D segments" << endl;
235 
236  if (nsegm!=0) {
237  // Find the best RecHit: look for the 4D RecHit with the phi angle closest
238  // to that of segment made of SimHits.
239  // RecHits must have delta alpha and delta position within 5 sigma of
240  // the residual distribution (we are looking for residuals of segments
241  // usefull to the track fit) for efficency purpose
242  const DTRecSegment4D* bestRecHit = 0;
243  bool bestRecHitFound = false;
244  double deltaAlpha = 99999;
245 
246  // Loop over the recHits of this chamberId
247  for (DTRecSegment4DCollection::const_iterator segment4D = range.first;
248  segment4D!=range.second;
249  ++segment4D){
250 
251  if (local) {
252  const DTChamberRecSegment2D* phiSeg = (*segment4D).phiSegment();
253  const DTSLRecSegment2D* zSeg = (*segment4D).zSegment();
254 
255  float t0phi = -999;
256  float t0z = -999;
257  int nHitPhi=0;
258  int nHitZ=0;
259  if (phiSeg) {
260  t0phi = phiSeg->t0();
261  nHitPhi = phiSeg->recHits().size();
262  }
263 
264  if (zSeg) {
265  t0z = zSeg->t0();
266  nHitZ = zSeg->recHits().size();
267  }
268 
269  hHitMult[abs((*chamberId).wheel())][(*chamberId).station()-1]->Fill(nHitPhi,nHitZ);
270  ht0[abs((*chamberId).wheel())][(*chamberId).station()-1]->Fill(t0phi,t0z);
271 
272  // Uncomment to skip segments w/o t0 computed
273  // if (!(t0phi>-998&&t0z>-998)) continue
274  }
275 
276  // Check the dimension
277  if((*segment4D).dimension() != 4) {
278  if(debug)cout << "[DTSegment4DQuality]***Error: This is not 4D segment!!!" << endl;
279  continue;
280  }
281  // Segment Local Direction and position (in Chamber RF)
282  LocalVector recSegDirection = (*segment4D).localDirection();
283  //LocalPoint recSegPosition = (*segment4D).localPosition();
284 
285  float recSegAlpha = DTHitQualityUtils::findSegmentAlphaAndBeta(recSegDirection).first;
286  if(debug)
287  cout << " RecSegment direction: " << recSegDirection << endl
288  << " position : " << (*segment4D).localPosition() << endl
289  << " alpha : " << recSegAlpha << endl
290  << " beta : " << DTHitQualityUtils::findSegmentAlphaAndBeta(recSegDirection).second << endl;
291 
292  if(fabs(recSegAlpha - alphaSimSeg) < deltaAlpha) {
293  deltaAlpha = fabs(recSegAlpha - alphaSimSeg);
294  bestRecHit = &(*segment4D);
295  bestRecHitFound = true;
296  }
297  } // End of Loop over all 4D RecHits
298 
299  if(bestRecHitFound) {
300  // Best rechit direction and position in Chamber RF
301  LocalPoint bestRecHitLocalPos = bestRecHit->localPosition();
302  LocalVector bestRecHitLocalDir = bestRecHit->localDirection();
303  // Errors on x and y
304  LocalError bestRecHitLocalPosErr = bestRecHit->localPositionError();
305  LocalError bestRecHitLocalDirErr = bestRecHit->localDirectionError();
306 
307  float alphaBestRHit = DTHitQualityUtils::findSegmentAlphaAndBeta(bestRecHitLocalDir).first;
308  float betaBestRHit = DTHitQualityUtils::findSegmentAlphaAndBeta(bestRecHitLocalDir).second;
309  // Errors on alpha and beta
310 
311  // Get position and direction using the rx projection (so in SL
312  // reference frame). Note that x (and y) are swapped wrt to Chamber
313  // frame
314  //if (bestRecHit->hasZed() ) {
315  const DTSLRecSegment2D * zedRecSeg = bestRecHit->zSegment();
316  LocalPoint bestRecHitLocalPosRZ = zedRecSeg->localPosition();
317  LocalVector bestRecHitLocalDirRZ = zedRecSeg->localDirection();
318  // Errors on x and y
319  LocalError bestRecHitLocalPosErrRZ = zedRecSeg->localPositionError();
320  LocalError bestRecHitLocalDirErrRZ = zedRecSeg->localDirectionError();
321 
322  float alphaBestRHitRZ = DTHitQualityUtils::findSegmentAlphaAndBeta(bestRecHitLocalDirRZ).first;
323 
324  // Get SimSeg position and Direction in rZ SL frame
325  const DTSuperLayer* sl = dtGeom->superLayer(zedRecSeg->superLayerId());
326  LocalPoint simSegLocalPosRZTmp = sl->toLocal(simSegmGlobalPos);
327  LocalVector simSegLocalDirRZ = sl->toLocal(simSegmGlobalDir);
328  LocalPoint simSegLocalPosRZ =
329  simSegLocalPosRZTmp + simSegLocalDirRZ*(-simSegLocalPosRZTmp.z()/(cos(simSegLocalDirRZ.theta())));
330  float alphaSimSegRZ = DTHitQualityUtils::findSegmentAlphaAndBeta(simSegLocalDirRZ).first;
331 
332  if (debug) cout <<
333  "RZ SL: recPos " << bestRecHitLocalPosRZ <<
334  "recDir " << bestRecHitLocalDirRZ <<
335  "recAlpha " << alphaBestRHitRZ << endl <<
336  "RZ SL: simPos " << simSegLocalPosRZ <<
337  "simDir " << simSegLocalDirRZ <<
338  "simAlpha " << alphaSimSegRZ << endl ;
339  //}
340 
341 
342  if(fabs(alphaBestRHit - alphaSimSeg) < 5*sigmaResAlpha &&
343  fabs(betaBestRHit - betaSimSeg) < 5*sigmaResBeta &&
344  fabs(bestRecHitLocalPos.x() - xSimSeg) < 5*sigmaResX &&
345  fabs(bestRecHitLocalPos.y() - ySimSeg) < 5*sigmaResY) {
346  recHitFound = true;
347  }
348 
349  // Fill Residual histos
350  HRes4DHit *histo=0;
351 
352  if((*chamberId).wheel() == 0)
353  histo = h4DHit_W0;
354  else if(abs((*chamberId).wheel()) == 1)
355  histo = h4DHit_W1;
356  else if(abs((*chamberId).wheel()) == 2)
357  histo = h4DHit_W2;
358 
359  histo->Fill(alphaSimSeg,
360  alphaBestRHit,
361  betaSimSeg,
362  betaBestRHit,
363  xSimSeg,
364  bestRecHitLocalPos.x(),
365  ySimSeg,
366  bestRecHitLocalPos.y(),
367  etaSimSeg,
368  phiSimSeg,
369  bestRecHitLocalPosRZ.x(),
370  simSegLocalPosRZ.x(),
371  alphaBestRHitRZ,
372  alphaSimSegRZ,
373  sqrt(bestRecHitLocalDirErr.xx()),
374  sqrt(bestRecHitLocalDirErr.yy()),
375  sqrt(bestRecHitLocalPosErr.xx()),
376  sqrt(bestRecHitLocalPosErr.yy()),
377  sqrt(bestRecHitLocalDirErrRZ.xx()),
378  sqrt(bestRecHitLocalPosErrRZ.xx())
379  );
380 
381  h4DHit->Fill(alphaSimSeg,
382  alphaBestRHit,
383  betaSimSeg,
384  betaBestRHit,
385  xSimSeg,
386  bestRecHitLocalPos.x(),
387  ySimSeg,
388  bestRecHitLocalPos.y(),
389  etaSimSeg,
390  phiSimSeg,
391  bestRecHitLocalPosRZ.x(),
392  simSegLocalPosRZ.x(),
393  alphaBestRHitRZ,
394  alphaSimSegRZ,
395  sqrt(bestRecHitLocalDirErr.xx()),
396  sqrt(bestRecHitLocalDirErr.yy()),
397  sqrt(bestRecHitLocalPosErr.xx()),
398  sqrt(bestRecHitLocalPosErr.yy()),
399  sqrt(bestRecHitLocalDirErrRZ.xx()),
400  sqrt(bestRecHitLocalPosErrRZ.xx())
401  );
402 
403  if (local) h4DHitWS[abs((*chamberId).wheel())][(*chamberId).station()-1]->Fill(alphaSimSeg,
404  alphaBestRHit,
405  betaSimSeg,
406  betaBestRHit,
407  xSimSeg,
408  bestRecHitLocalPos.x(),
409  ySimSeg,
410  bestRecHitLocalPos.y(),
411  etaSimSeg,
412  phiSimSeg,
413  bestRecHitLocalPosRZ.x(),
414  simSegLocalPosRZ.x(),
415  alphaBestRHitRZ,
416  alphaSimSegRZ,
417  sqrt(bestRecHitLocalDirErr.xx()),
418  sqrt(bestRecHitLocalDirErr.yy()),
419  sqrt(bestRecHitLocalPosErr.xx()),
420  sqrt(bestRecHitLocalPosErr.yy()),
421  sqrt(bestRecHitLocalDirErrRZ.xx()),
422  sqrt(bestRecHitLocalPosErrRZ.xx())
423  );
424 
425 
426  } //end of if(bestRecHitFound)
427 
428  } //end of if(nsegm!=0)
429 
430  // Fill Efficiency plot
431  if(doall){
432  HEff4DHit *heff = 0;
433 
434  if((*chamberId).wheel() == 0)
435  heff = hEff_W0;
436  else if(abs((*chamberId).wheel()) == 1)
437  heff = hEff_W1;
438  else if(abs((*chamberId).wheel()) == 2)
439  heff = hEff_W2;
440  heff->Fill(etaSimSeg, phiSimSeg, xSimSeg, ySimSeg, alphaSimSeg, betaSimSeg, recHitFound);
441  hEff_All->Fill(etaSimSeg, phiSimSeg, xSimSeg, ySimSeg, alphaSimSeg, betaSimSeg, recHitFound);
442  if (local) hEffWS[abs((*chamberId).wheel())][(*chamberId).station()-1]->Fill(etaSimSeg, phiSimSeg, xSimSeg, ySimSeg, alphaSimSeg, betaSimSeg, recHitFound);
443 
444  }
445  } // End of loop over chambers
446  }
virtual LocalError localPositionError() const
local position error in SL frame
float xx() const
Definition: LocalError.h:24
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:52
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:47
Geom::Phi< T > phi() const
Definition: PV3DBase.h:68
T y() const
Definition: PV3DBase.h:62
#define abs(x)
Definition: mlp_lapack.h:159
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:62
identifier iterator
Definition: RangeMap.h:138
virtual LocalError localDirectionError() const
the local direction error (xx,xy,yy) in SL frame: only xx is not 0.
static std::map< DTWireId, const PSimHit * > mapMuSimHitsPerWire(const std::map< DTWireId, edm::PSimHitContainer > &simHitWireMap)
Create a map between the Mu SimHits and corresponding MuBarWireId ;.
void Fill(long long x)
Geom::Theta< T > theta() const
Definition: PV3DBase.h:74
virtual LocalError localPositionError() const
Local position error in Chamber frame.
static std::pair< const PSimHit *, const PSimHit * > findMuSimSegment(const std::map< DTWireId, const PSimHit * > &mapWireAndMuSimHit)
Find Innermost and outermost SimHit from Mu in a SL (they identify a simulated segment) ...
virtual LocalVector localDirection() const
Local direction in Chamber frame.
edm::InputTag segment4DLabel
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:45
HRes4DHit * h4DHitWS[3][4]
float yy() const
Definition: LocalError.h:26
virtual std::vector< const TrackingRecHit * > recHits() const
Access to component RecHits (if any)
T sqrt(T t)
Definition: SSEVec.h:46
void Fill(float simDirectionAlpha, float recDirectionAlpha, float simDirectionBeta, float recDirectionBeta, float simX, float recX, float simY, float recY, float simEta, float simPhi, float recYRZ, float simYRZ, float recBetaRZ, float simBetaRZ, float sigmaAlpha, float sigmaBeta, float sigmaX, float sigmaY, float sigmaBetaRZ, float sigmaYRZ)
Definition: Histograms.h:903
T z() const
Definition: PV3DBase.h:63
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
MonitorElement * hHitMult[3][4]
A set of histograms for efficiency 4D RecHits.
Definition: Histograms.h:1116
bool isValid() const
Definition: HandleBase.h:76
virtual LocalPoint localPosition() const
Local position in Chamber frame.
static std::pair< double, double > findSegmentAlphaAndBeta(const LocalVector &direction)
Find the angles from a segment direction:
DTSuperLayerId superLayerId() const
The id of the superlayer on which reside the segment.
virtual LocalPoint localPosition() const
local position in SL frame
const DTSLRecSegment2D * zSegment() const
The Z segment: 0 if not zed projection available.
static std::pair< LocalVector, LocalPoint > findMuSimSegmentDirAndPos(const std::pair< const PSimHit *, const PSimHit * > &inAndOutSimHit, const DetId detId, const DTGeometry *muonGeom)
Find direction and position of a segment (in local RF) from outer and inner mu SimHit in the RF of ob...
edm::InputTag simHitLabel
const T & get() const
Definition: EventSetup.h:55
tuple simHits
Definition: trackerHits.py:16
HEff4DHit * hEffWS[3][4]
T eta() const
Definition: PV3DBase.h:75
virtual LocalVector localDirection() const
the local direction in SL frame
tuple cout
Definition: gather_cfg.py:121
std::vector< PSimHit > PSimHitContainer
static std::map< DTWireId, edm::PSimHitContainer > mapSimHitsPerWire(const edm::PSimHitContainer &simhits)
double t0() const
Get the segment t0 (if recomputed, 0 is returned otherwise)
T x() const
Definition: PV3DBase.h:61
void Fill(float etaSimSegm, float phiSimSegm, float xSimSegm, float ySimSegm, float alphaSimSegm, float betaSimSegm, bool fillRecHit)
Definition: Histograms.h:1211
MonitorElement * ht0[3][4]
virtual LocalError localDirectionError() const
Local direction error in the Chamber frame.
void DTSegment4DQuality::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 117 of file DTSegment4DQuality.cc.

117  {
118  // Write the histos to file
119  //theFile->cd();
120 
121  //h4DHit->Write();
122  //h4DHit_W0->Write();
123  //h4DHit_W1->Write();
124  //h4DHit_W2->Write();
125 
126  if(doall){
131  }
132  //hEff_All->Write();
133  //hEff_W0->Write();
134  //hEff_W1->Write();
135  //hEff_W2->Write();
136  //if ( rootFileName.size() != 0 && dbe_ ) dbe_->save(rootFileName);
137 
138  //theFile->Close();
139 }
void ComputeEfficiency()
Definition: Histograms.h:1238
void DTSegment4DQuality::endLuminosityBlock ( edm::LuminosityBlock const &  lumiSeg,
edm::EventSetup const &  c 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 111 of file DTSegment4DQuality.cc.

112  {
113 
114 
115 }

Member Data Documentation

DQMStore* DTSegment4DQuality::dbe_
private

Definition at line 95 of file DTSegment4DQuality.h.

bool DTSegment4DQuality::debug
private

Definition at line 67 of file DTSegment4DQuality.h.

bool DTSegment4DQuality::doall
private

Definition at line 96 of file DTSegment4DQuality.h.

HRes4DHit* DTSegment4DQuality::h4DHit
private

Definition at line 80 of file DTSegment4DQuality.h.

HRes4DHit* DTSegment4DQuality::h4DHit_W0
private

Definition at line 81 of file DTSegment4DQuality.h.

HRes4DHit* DTSegment4DQuality::h4DHit_W1
private

Definition at line 82 of file DTSegment4DQuality.h.

HRes4DHit* DTSegment4DQuality::h4DHit_W2
private

Definition at line 83 of file DTSegment4DQuality.h.

HRes4DHit* DTSegment4DQuality::h4DHitWS[3][4]
private

Definition at line 84 of file DTSegment4DQuality.h.

HEff4DHit* DTSegment4DQuality::hEff_All
private

Definition at line 86 of file DTSegment4DQuality.h.

HEff4DHit* DTSegment4DQuality::hEff_W0
private

Definition at line 87 of file DTSegment4DQuality.h.

HEff4DHit* DTSegment4DQuality::hEff_W1
private

Definition at line 88 of file DTSegment4DQuality.h.

HEff4DHit* DTSegment4DQuality::hEff_W2
private

Definition at line 89 of file DTSegment4DQuality.h.

HEff4DHit* DTSegment4DQuality::hEffWS[3][4]
private

Definition at line 90 of file DTSegment4DQuality.h.

MonitorElement* DTSegment4DQuality::hHitMult[3][4]
private

Definition at line 92 of file DTSegment4DQuality.h.

MonitorElement* DTSegment4DQuality::ht0[3][4]
private

Definition at line 93 of file DTSegment4DQuality.h.

bool DTSegment4DQuality::local
private

Definition at line 97 of file DTSegment4DQuality.h.

std::string DTSegment4DQuality::rootFileName
private

Definition at line 69 of file DTSegment4DQuality.h.

edm::InputTag DTSegment4DQuality::segment4DLabel
private

Definition at line 72 of file DTSegment4DQuality.h.

double DTSegment4DQuality::sigmaResAlpha
private

Definition at line 77 of file DTSegment4DQuality.h.

double DTSegment4DQuality::sigmaResBeta
private

Definition at line 78 of file DTSegment4DQuality.h.

double DTSegment4DQuality::sigmaResX
private

Definition at line 74 of file DTSegment4DQuality.h.

double DTSegment4DQuality::sigmaResY
private

Definition at line 75 of file DTSegment4DQuality.h.

edm::InputTag DTSegment4DQuality::simHitLabel
private

Definition at line 71 of file DTSegment4DQuality.h.