CMS 3D CMS Logo

TrackSplittingMonitor.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author Suchandra Dutta , Giorgia Mila
5  */
6 
15 
19 //#include "DQM/TrackingMonitor/interface/TrackAnalyzer.h"
21 
24 
30 #include <string>
31 
33  : dqmStore_( edm::Service<DQMStore>().operator->() )
34  , conf_( iConfig )
35 {
36 
37  splitTracks_ = conf_.getParameter< edm::InputTag >("splitTrackCollection");
38  splitMuons_ = conf_.getParameter< edm::InputTag >("splitMuonCollection");
39  splitTracksToken_ = consumes<std::vector<reco::Track> >(splitTracks_);
40  splitMuonsToken_ = mayConsume<std::vector<reco::Muon> >(splitMuons_);
41 
42  plotMuons_ = conf_.getParameter<bool>("ifPlotMuons");
43 
44  // cuts
45  pixelHitsPerLeg_ = conf_.getParameter<int>("pixelHitsPerLeg");
46  totalHitsPerLeg_ = conf_.getParameter<int>("totalHitsPerLeg");
47  d0Cut_ = conf_.getParameter<double>("d0Cut");
48  dzCut_ = conf_.getParameter<double>("dzCut");
49  ptCut_ = conf_.getParameter<double>("ptCut");
50  norchiCut_ = conf_.getParameter<double>("norchiCut");
51 
52 }
53 
55 }
56 
58  edm::Run const & /* iRun */,
59  edm::EventSetup const & /* iSetup */)
60 
61 {
62 
63  std::string MEFolderName = conf_.getParameter<std::string>("FolderName");
64  ibooker.setCurrentFolder(MEFolderName);
65 
66  // bin declarations
67  int ddxyBin = conf_.getParameter<int>("ddxyBin");
68  double ddxyMin = conf_.getParameter<double>("ddxyMin");
69  double ddxyMax = conf_.getParameter<double>("ddxyMax");
70 
71  int ddzBin = conf_.getParameter<int>("ddzBin");
72  double ddzMin = conf_.getParameter<double>("ddzMin");
73  double ddzMax = conf_.getParameter<double>("ddzMax");
74 
75  int dphiBin = conf_.getParameter<int>("dphiBin");
76  double dphiMin = conf_.getParameter<double>("dphiMin");
77  double dphiMax = conf_.getParameter<double>("dphiMax");
78 
79  int dthetaBin = conf_.getParameter<int>("dthetaBin");
80  double dthetaMin = conf_.getParameter<double>("dthetaMin");
81  double dthetaMax = conf_.getParameter<double>("dthetaMax");
82 
83  int dptBin = conf_.getParameter<int>("dptBin");
84  double dptMin = conf_.getParameter<double>("dptMin");
85  double dptMax = conf_.getParameter<double>("dptMax");
86 
87  int dcurvBin = conf_.getParameter<int>("dcurvBin");
88  double dcurvMin = conf_.getParameter<double>("dcurvMin");
89  double dcurvMax = conf_.getParameter<double>("dcurvMax");
90 
91  int normBin = conf_.getParameter<int>("normBin");
92  double normMin = conf_.getParameter<double>("normMin");
93  double normMax = conf_.getParameter<double>("normMax");
94 
95  // declare histogram
96  ddxyAbsoluteResiduals_tracker_ = ibooker.book1D( "ddxyAbsoluteResiduals_tracker", "ddxyAbsoluteResiduals_tracker", ddxyBin, ddxyMin, ddxyMax );
97  ddzAbsoluteResiduals_tracker_ = ibooker.book1D( "ddzAbsoluteResiduals_tracker", "ddzAbsoluteResiduals_tracker", ddzBin, ddzMin, ddzMax );
98  dphiAbsoluteResiduals_tracker_ = ibooker.book1D( "dphiAbsoluteResiduals_tracker", "dphiAbsoluteResiduals_tracker", dphiBin, dphiMin, dphiMax );
99  dthetaAbsoluteResiduals_tracker_ = ibooker.book1D( "dthetaAbsoluteResiduals_tracker", "dthetaAbsoluteResiduals_tracker", dthetaBin, dthetaMin, dthetaMax );
100  dptAbsoluteResiduals_tracker_ = ibooker.book1D( "dptAbsoluteResiduals_tracker", "dptAbsoluteResiduals_tracker", dptBin, dptMin, dptMax );
101  dcurvAbsoluteResiduals_tracker_ = ibooker.book1D( "dcurvAbsoluteResiduals_tracker", "dcurvAbsoluteResiduals_tracker", dcurvBin, dcurvMin, dcurvMax );
102 
103  ddxyNormalizedResiduals_tracker_ = ibooker.book1D( "ddxyNormalizedResiduals_tracker", "ddxyNormalizedResiduals_tracker", normBin, normMin, normMax );
104  ddzNormalizedResiduals_tracker_ = ibooker.book1D( "ddzNormalizedResiduals_tracker", "ddzNormalizedResiduals_tracker", normBin, normMin, normMax );
105  dphiNormalizedResiduals_tracker_ = ibooker.book1D( "dphiNormalizedResiduals_tracker", "dphiNormalizedResiduals_tracker", normBin, normMin, normMax );
106  dthetaNormalizedResiduals_tracker_ = ibooker.book1D( "dthetaNormalizedResiduals_tracker", "dthetaNormalizedResiduals_tracker", normBin, normMin, normMax );
107  dptNormalizedResiduals_tracker_ = ibooker.book1D( "dptNormalizedResiduals_tracker", "dptNormalizedResiduals_tracker", normBin, normMin, normMax );
108  dcurvNormalizedResiduals_tracker_ = ibooker.book1D( "dcurvNormalizedResiduals_tracker", "dcurvNormalizedResiduals_tracker", normBin, normMin, normMax );
109 
110  if (plotMuons_){
111  ddxyAbsoluteResiduals_global_ = ibooker.book1D( "ddxyAbsoluteResiduals_global", "ddxyAbsoluteResiduals_global", ddxyBin, ddxyMin, ddxyMax );
112  ddzAbsoluteResiduals_global_ = ibooker.book1D( "ddzAbsoluteResiduals_global", "ddzAbsoluteResiduals_global", ddzBin, ddzMin, ddzMax );
113  dphiAbsoluteResiduals_global_ = ibooker.book1D( "dphiAbsoluteResiduals_global", "dphiAbsoluteResiduals_global", dphiBin, dphiMin, dphiMax );
114  dthetaAbsoluteResiduals_global_ = ibooker.book1D( "dthetaAbsoluteResiduals_global", "dthetaAbsoluteResiduals_global", dthetaBin, dthetaMin, dthetaMax );
115  dptAbsoluteResiduals_global_ = ibooker.book1D( "dptAbsoluteResiduals_global", "dptAbsoluteResiduals_global", dptBin, dptMin, dptMax );
116  dcurvAbsoluteResiduals_global_ = ibooker.book1D( "dcurvAbsoluteResiduals_global", "dcurvAbsoluteResiduals_global", dcurvBin, dcurvMin, dcurvMax );
117 
118  ddxyNormalizedResiduals_global_ = ibooker.book1D( "ddxyNormalizedResiduals_global", "ddxyNormalizedResiduals_global", normBin, normMin, normMax );
119  ddzNormalizedResiduals_global_ = ibooker.book1D( "ddzNormalizedResiduals_global", "ddzNormalizedResiduals_global", normBin, normMin, normMax );
120  dphiNormalizedResiduals_global_ = ibooker.book1D( "dphiNormalizedResiduals_global", "dphiNormalizedResiduals_global", normBin, normMin, normMax );
121  dthetaNormalizedResiduals_global_ = ibooker.book1D( "dthetaNormalizedResiduals_global", "dthetaNormalizedResiduals_global", normBin, normMin, normMax );
122  dptNormalizedResiduals_global_ = ibooker.book1D( "dptNormalizedResiduals_global", "dptNormalizedResiduals_global", normBin, normMin, normMax );
123  dcurvNormalizedResiduals_global_ = ibooker.book1D( "dcurvNormalizedResiduals_global", "dcurvNormalizedResiduals_global", normBin, normMin, normMax );
124  }
125 
126  ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta d_{xy})/#sqrt{2} [#mum]" );
127  ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta d_{z})/#sqrt{2} [#mum]" );
128  ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta #phi)/#sqrt{2} [mrad]" );
129  ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta #theta)/#sqrt{2} [mrad]" );
130  ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta pT)/#sqrt{2} [GeV]" );
131  ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta (1/pT))/#sqrt{2} [GeV^{-1}]" );
132 
133  ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta d_{xy}/#sigma(d_{xy}" );
134  ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta d_{z}/#sigma(d_{z})" );
135  ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta #phi/#sigma(d_{#phi})" );
136  ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta #theta/#sigma(d_{#theta})" );
137  ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta p_{T}/#sigma(p_{T})" );
138  ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta 1/p_{T}/#sigma(1/p_{T})" );
139 
140  if (plotMuons_){
141  ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta d_{xy})/#sqrt{2} [#mum]" );
142  ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta d_{z})/#sqrt{2} [#mum]" );
143  ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta #phi)/#sqrt{2} [mrad]" );
144  ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta #theta)/#sqrt{2} [mrad]" );
145  ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta pT)/#sqrt{2} [GeV]" );
146  ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta (1/pT))/#sqrt{2} [GeV^{-1}]" );
147 
148  ddxyNormalizedResiduals_global_->setAxisTitle( "#delta d_{xy}/#sigma(d_{xy}" );
149  ddxyNormalizedResiduals_global_->setAxisTitle( "#delta d_{z}/#sigma(d_{z})" );
150  ddxyNormalizedResiduals_global_->setAxisTitle( "#delta #phi/#sigma(d_{#phi})" );
151  ddxyNormalizedResiduals_global_->setAxisTitle( "#delta #theta/#sigma(d_{#theta})" );
152  ddxyNormalizedResiduals_global_->setAxisTitle( "#delta p_{T}/#sigma(p_{T})" );
153  ddxyNormalizedResiduals_global_->setAxisTitle( "#delta 1/p_{T}/#sigma(1/p_{T})" );
154  }
155 
156 }
157 
158 
159 //
160 // -- Analyse
161 //
163 
164 
165  iSetup.get<IdealMagneticFieldRecord>().get(theMagField);
167  iSetup.get<MuonGeometryRecord>().get(dtGeometry);
168  iSetup.get<MuonGeometryRecord>().get(cscGeometry);
169  iSetup.get<MuonGeometryRecord>().get(rpcGeometry);
170 
172  iEvent.getByToken(splitTracksToken_, splitTracks);
173  if (!splitTracks.isValid()) return;
174 
176  if (plotMuons_){
177  iEvent.getByToken(splitMuonsToken_, splitMuons);
178  }
179 
180  if (splitTracks->size() == 2){
181  // check that there are 2 tracks in split track collection
182  edm::LogInfo("TrackSplittingMonitor") << "Split Track size: " << splitTracks->size();
183 
184  // split tracks calculations
185  reco::Track track1 = splitTracks->at(0);
186  reco::Track track2 = splitTracks->at(1);
187 
188 
189  // -------------------------- basic selection ---------------------------
190 
191  // hit counting
192  // looping through the hits for track 1
193  double nRechits1 =0;
194  double nRechitinBPIX1 =0;
195  for (trackingRecHit_iterator iHit = track1.recHitsBegin(); iHit != track1.recHitsEnd(); ++iHit) {
196  if((*iHit)->isValid()) {
197  nRechits1++;
198  int type =(*iHit)->geographicalId().subdetId();
199  if(type==int(PixelSubdetector::PixelBarrel)){++nRechitinBPIX1;}
200  }
201  }
202  // looping through the hits for track 2
203  double nRechits2 =0;
204  double nRechitinBPIX2 =0;
205  for (trackingRecHit_iterator iHit = track2.recHitsBegin(); iHit != track2.recHitsEnd(); ++iHit) {
206  if((*iHit)->isValid()) {
207  nRechits2++;
208  int type =(*iHit)->geographicalId().subdetId();
209  if(type==int(PixelSubdetector::PixelBarrel)){++nRechitinBPIX2;}
210  }
211  }
212 
213  // DCA of each track
214  double d01 = track1.d0();
215  double dz1 = track1.dz();
216  double d02 = track2.d0();
217  double dz2 = track2.dz();
218 
219  // pT of each track
220  double pt1 = track1.pt();
221  double pt2 = track2.pt();
222 
223  // chi2 of each track
224  double norchi1 = track1.normalizedChi2();
225  double norchi2 = track2.normalizedChi2();
226 
227  // basic selection
228  // pixel hits and total hits
229  if ((nRechitinBPIX1 >= pixelHitsPerLeg_)&&(nRechitinBPIX1 >= pixelHitsPerLeg_)&&(nRechits1 >= totalHitsPerLeg_)&&(nRechits2 >= totalHitsPerLeg_)){
230  // dca cut
231  if ( ((fabs(d01) < d0Cut_))&&(fabs(d02) < d0Cut_)&&(fabs(dz2) < dzCut_)&&(fabs(dz2) < dzCut_) ){
232  // pt cut
233  if ( (pt1+pt2)/2 < ptCut_){
234  // chi2 cut
235  if ( (norchi1 < norchiCut_)&&(norchi2 < norchiCut_) ){
236 
237  // passed all cuts...
238  edm::LogInfo("TrackSplittingMonitor") << " Setected after all cuts ?";
239 
240  double ddxyVal = d01 - d02;
241  double ddzVal = dz1 - dz2;
242  double dphiVal = track1.phi() - track2.phi();
243  double dthetaVal = track1.theta() - track2.theta();
244  double dptVal = pt1 - pt2;
245  double dcurvVal = (1/pt1) - (1/pt2);
246 
247  double d01ErrVal = track1.d0Error();
248  double d02ErrVal = track2.d0Error();
249  double dz1ErrVal = track1.dzError();
250  double dz2ErrVal = track2.dzError();
251  double phi1ErrVal = track1.phiError();
252  double phi2ErrVal = track2.phiError();
253  double theta1ErrVal = track1.thetaError();
254  double theta2ErrVal = track2.thetaError();
255  double pt1ErrVal = track1.ptError();
256  double pt2ErrVal = track2.ptError();
257 
258  ddxyAbsoluteResiduals_tracker_->Fill( 10000.0*ddxyVal/sqrt(2.0) );
259  ddxyAbsoluteResiduals_tracker_->Fill( 10000.0*ddzVal/sqrt(2.0) );
260  ddxyAbsoluteResiduals_tracker_->Fill( 1000.0*dphiVal/sqrt(2.0) );
261  ddxyAbsoluteResiduals_tracker_->Fill( 1000.0*dthetaVal/sqrt(2.0) );
262  ddxyAbsoluteResiduals_tracker_->Fill( dptVal/sqrt(2.0) );
263  ddxyAbsoluteResiduals_tracker_->Fill( dcurvVal/sqrt(2.0) );
264 
265  ddxyNormalizedResiduals_tracker_->Fill( ddxyVal/sqrt( d01ErrVal*d01ErrVal + d02ErrVal*d02ErrVal ) );
266  ddxyNormalizedResiduals_tracker_->Fill( ddzVal/sqrt( dz1ErrVal*dz1ErrVal + dz2ErrVal*dz2ErrVal ) );
267  ddxyNormalizedResiduals_tracker_->Fill( dphiVal/sqrt( phi1ErrVal*phi1ErrVal + phi2ErrVal*phi2ErrVal ) );
268  ddxyNormalizedResiduals_tracker_->Fill( dthetaVal/sqrt( theta1ErrVal*theta1ErrVal + theta2ErrVal*theta2ErrVal ) );
269  ddxyNormalizedResiduals_tracker_->Fill( dptVal/sqrt( pt1ErrVal*pt1ErrVal + pt2ErrVal*pt2ErrVal ) );
270  ddxyNormalizedResiduals_tracker_->Fill( dcurvVal/sqrt( pow(pt1ErrVal,2)/pow(pt1,4) + pow(pt2ErrVal,2)/pow(pt2,4) ) );
271 
272  // if do the same for split muons
273  if (plotMuons_ && splitMuons.isValid()){
274 
275  int gmCtr = 0;
276  bool topGlobalMuonFlag = false;
277  bool bottomGlobalMuonFlag = false;
278  int topGlobalMuon = -1;
279  int bottomGlobalMuon = -1;
280  double topGlobalMuonNorchi2 = 1e10;
281  double bottomGlobalMuonNorchi2 = 1e10;
282 
283  // check if usable split global muons
284  for (std::vector<reco::Muon>::const_iterator gmI = splitMuons->begin(); gmI != splitMuons->end(); gmI++){
285  if ( gmI->isTrackerMuon() && gmI->isStandAloneMuon() && gmI->isGlobalMuon() ){
286 
287  reco::TrackRef trackerTrackRef1( splitTracks, 0 );
288  reco::TrackRef trackerTrackRef2( splitTracks, 1 );
289 
290  if (gmI->innerTrack() == trackerTrackRef1){
291  if (gmI->globalTrack()->normalizedChi2() < topGlobalMuonNorchi2){
292  topGlobalMuonFlag = true;
293  topGlobalMuonNorchi2 = gmI->globalTrack()->normalizedChi2();
294  topGlobalMuon = gmCtr;
295  }
296  }
297  if (gmI->innerTrack() == trackerTrackRef2){
298  if (gmI->globalTrack()->normalizedChi2() < bottomGlobalMuonNorchi2){
299  bottomGlobalMuonFlag = true;
300  bottomGlobalMuonNorchi2 = gmI->globalTrack()->normalizedChi2();
301  bottomGlobalMuon = gmCtr;
302  }
303  }
304  }
305  gmCtr++;
306  }
307 
308  if (bottomGlobalMuonFlag && topGlobalMuonFlag) {
309 
310  reco::Muon muonTop = splitMuons->at( topGlobalMuon );
311  reco::Muon muonBottom = splitMuons->at( bottomGlobalMuon );
312 
313  reco::TrackRef glb1 = muonTop.globalTrack();
314  reco::TrackRef glb2 = muonBottom.globalTrack();
315 
316  double ddxyValGlb = glb1->d0() - glb2->d0();
317  double ddzValGlb = glb1->dz() - glb2->dz();
318  double dphiValGlb = glb1->phi() - glb2->phi();
319  double dthetaValGlb = glb1->theta() - glb2->theta();
320  double dptValGlb = glb1->pt() - glb2->pt();
321  double dcurvValGlb = (1/glb1->pt()) - (1/glb2->pt());
322 
323  double d01ErrValGlb = glb1->d0Error();
324  double d02ErrValGlb = glb2->d0Error();
325  double dz1ErrValGlb = glb1->dzError();
326  double dz2ErrValGlb = glb2->dzError();
327  double phi1ErrValGlb = glb1->phiError();
328  double phi2ErrValGlb = glb2->phiError();
329  double theta1ErrValGlb = glb1->thetaError();
330  double theta2ErrValGlb = glb2->thetaError();
331  double pt1ErrValGlb = glb1->ptError();
332  double pt2ErrValGlb = glb2->ptError();
333 
334  ddxyAbsoluteResiduals_global_->Fill( 10000.0*ddxyValGlb/sqrt(2.0) );
335  ddxyAbsoluteResiduals_global_->Fill( 10000.0*ddzValGlb/sqrt(2.0) );
336  ddxyAbsoluteResiduals_global_->Fill( 1000.0*dphiValGlb/sqrt(2.0) );
337  ddxyAbsoluteResiduals_global_->Fill( 1000.0*dthetaValGlb/sqrt(2.0) );
338  ddxyAbsoluteResiduals_global_->Fill( dptValGlb/sqrt(2.0) );
339  ddxyAbsoluteResiduals_global_->Fill( dcurvValGlb/sqrt(2.0) );
340 
341  ddxyNormalizedResiduals_global_->Fill( ddxyValGlb/sqrt( d01ErrValGlb*d01ErrValGlb + d02ErrValGlb*d02ErrValGlb ) );
342  ddxyNormalizedResiduals_global_->Fill( ddzValGlb/sqrt( dz1ErrValGlb*dz1ErrValGlb + dz2ErrValGlb*dz2ErrValGlb ) );
343  ddxyNormalizedResiduals_global_->Fill( dphiValGlb/sqrt( phi1ErrValGlb*phi1ErrValGlb + phi2ErrValGlb*phi2ErrValGlb ) );
344  ddxyNormalizedResiduals_global_->Fill( dthetaValGlb/sqrt( theta1ErrValGlb*theta1ErrValGlb + theta2ErrValGlb*theta2ErrValGlb ) );
345  ddxyNormalizedResiduals_global_->Fill( dptValGlb/sqrt( pt1ErrValGlb*pt1ErrValGlb + pt2ErrValGlb*pt2ErrValGlb ) );
346  ddxyNormalizedResiduals_global_->Fill( dcurvValGlb/sqrt( pow(pt1ErrValGlb,2)/pow(pt1,4) + pow(pt2ErrValGlb,2)/pow(pt2,4) ) );
347 
348  }
349 
350 
351  } // end of split muons loop
352  }
353  }
354  }
355  }
356  }
357 }
358 
359 
360 
type
Definition: HCALResponse.h:21
edm::ESHandle< TrackerGeometry > theGeometry
T getParameter(std::string const &) const
double d0Error() const
error on d0
Definition: TrackBase.h:797
MonitorElement * dthetaAbsoluteResiduals_global_
MonitorElement * dcurvNormalizedResiduals_global_
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:592
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:556
MonitorElement * dthetaAbsoluteResiduals_tracker_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
double theta() const
polar angle
Definition: TrackBase.h:574
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
MonitorElement * ddzNormalizedResiduals_global_
MonitorElement * ddxyNormalizedResiduals_tracker_
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:640
MonitorElement * dcurvAbsoluteResiduals_tracker_
edm::ESHandle< RPCGeometry > rpcGeometry
MonitorElement * ddzAbsoluteResiduals_global_
edm::EDGetTokenT< std::vector< reco::Muon > > splitMuonsToken_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
void Fill(long long x)
MonitorElement * dphiNormalizedResiduals_global_
MonitorElement * dthetaNormalizedResiduals_tracker_
int iEvent
Definition: GenABIO.cc:230
MonitorElement * ddxyNormalizedResiduals_global_
MonitorElement * dthetaNormalizedResiduals_global_
MonitorElement * dphiNormalizedResiduals_tracker_
T sqrt(T t)
Definition: SSEVec.h:18
double pt() const
track transverse momentum
Definition: TrackBase.h:616
MonitorElement * dcurvAbsoluteResiduals_global_
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:758
MonitorElement * dphiAbsoluteResiduals_tracker_
double phiError() const
error on phi
Definition: TrackBase.h:785
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
edm::ESHandle< MagneticField > theMagField
MonitorElement * ddzNormalizedResiduals_tracker_
MonitorElement * dptAbsoluteResiduals_global_
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:104
bool isValid() const
Definition: HandleBase.h:74
MonitorElement * dptAbsoluteResiduals_tracker_
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:604
double dzError() const
error on dz
Definition: TrackBase.h:809
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
const T & get() const
Definition: EventSetup.h:56
TrackSplittingMonitor(const edm::ParameterSet &)
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
MonitorElement * ddzAbsoluteResiduals_tracker_
MonitorElement * ddxyAbsoluteResiduals_global_
HLT enums.
edm::ESHandle< DTGeometry > dtGeometry
MonitorElement * dcurvNormalizedResiduals_tracker_
MonitorElement * dphiAbsoluteResiduals_global_
edm::EDGetTokenT< std::vector< reco::Track > > splitTracksToken_
MonitorElement * dptNormalizedResiduals_tracker_
MonitorElement * dptNormalizedResiduals_global_
edm::ESHandle< CSCGeometry > cscGeometry
MonitorElement * ddxyAbsoluteResiduals_tracker_
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
Definition: Run.h:42
double thetaError() const
error on theta
Definition: TrackBase.h:767
virtual TrackRef globalTrack() const
reference to Track reconstructed in both tracked and muon detector
Definition: Muon.h:54
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:109