CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackBuildingAnalyzer.cc
Go to the documentation of this file.
11 
16 #include <string>
17 #include "TMath.h"
18 
19 #include <iostream>
20 
22  : conf_( iConfig )
23  , SeedPt(NULL)
24  , SeedEta(NULL)
25  , SeedPhi(NULL)
26  , SeedPhiVsEta(NULL)
27  , SeedTheta(NULL)
28  , SeedQ(NULL)
29  , SeedDxy(NULL)
30  , SeedDz(NULL)
31  , NumberOfRecHitsPerSeed(NULL)
32  , NumberOfRecHitsPerSeedVsPhiProfile(NULL)
33  , NumberOfRecHitsPerSeedVsEtaProfile(NULL)
34 {
35 }
36 
38 {
39 }
40 
42 {
43 
44  // parameters from the configuration
46  std::string MEFolderName = conf_.getParameter<std::string>("FolderName");
47 
48  // std::cout << "[TrackBuildingAnalyzer::beginRun] AlgoName: " << AlgoName << std::endl;
49 
50  // use the AlgoName and Quality Name
51  std::string CatagoryName = AlgoName;
52 
53  // get binning from the configuration
54  int TrackPtBin = conf_.getParameter<int>( "TrackPtBin");
55  double TrackPtMin = conf_.getParameter<double>("TrackPtMin");
56  double TrackPtMax = conf_.getParameter<double>("TrackPtMax");
57 
58  int PhiBin = conf_.getParameter<int>( "PhiBin");
59  double PhiMin = conf_.getParameter<double>("PhiMin");
60  double PhiMax = conf_.getParameter<double>("PhiMax");
61 
62  int EtaBin = conf_.getParameter<int>( "EtaBin");
63  double EtaMin = conf_.getParameter<double>("EtaMin");
64  double EtaMax = conf_.getParameter<double>("EtaMax");
65 
66  int ThetaBin = conf_.getParameter<int>( "ThetaBin");
67  double ThetaMin = conf_.getParameter<double>("ThetaMin");
68  double ThetaMax = conf_.getParameter<double>("ThetaMax");
69 
70  int TrackQBin = conf_.getParameter<int>( "TrackQBin");
71  double TrackQMin = conf_.getParameter<double>("TrackQMin");
72  double TrackQMax = conf_.getParameter<double>("TrackQMax");
73 
74  int SeedDxyBin = conf_.getParameter<int>( "SeedDxyBin");
75  double SeedDxyMin = conf_.getParameter<double>("SeedDxyMin");
76  double SeedDxyMax = conf_.getParameter<double>("SeedDxyMax");
77 
78  int SeedDzBin = conf_.getParameter<int>( "SeedDzBin");
79  double SeedDzMin = conf_.getParameter<double>("SeedDzMin");
80  double SeedDzMax = conf_.getParameter<double>("SeedDzMax");
81 
82  int SeedHitBin = conf_.getParameter<int>( "SeedHitBin");
83  double SeedHitMin = conf_.getParameter<double>("SeedHitMin");
84  double SeedHitMax = conf_.getParameter<double>("SeedHitMax");
85 
86  int TCDxyBin = conf_.getParameter<int>( "TCDxyBin");
87  double TCDxyMin = conf_.getParameter<double>("TCDxyMin");
88  double TCDxyMax = conf_.getParameter<double>("TCDxyMax");
89 
90  int TCDzBin = conf_.getParameter<int>( "TCDzBin");
91  double TCDzMin = conf_.getParameter<double>("TCDzMin");
92  double TCDzMax = conf_.getParameter<double>("TCDzMax");
93 
94  int TCHitBin = conf_.getParameter<int>( "TCHitBin");
95  double TCHitMin = conf_.getParameter<double>("TCHitMin");
96  double TCHitMax = conf_.getParameter<double>("TCHitMax");
97 
98 
99  edm::InputTag seedProducer = conf_.getParameter<edm::InputTag>("SeedProducer");
100  edm::InputTag tcProducer = conf_.getParameter<edm::InputTag>("TCProducer");
101 
102  doAllPlots = conf_.getParameter<bool>("doAllPlots");
103  doAllSeedPlots = conf_.getParameter<bool>("doSeedParameterHistos");
104  doTCPlots = conf_.getParameter<bool>("doTrackCandHistos");
105  doAllTCPlots = conf_.getParameter<bool>("doAllTrackCandHistos");
106  doPT = conf_.getParameter<bool>("doSeedPTHisto");
107  doETA = conf_.getParameter<bool>("doSeedETAHisto");
108  doPHI = conf_.getParameter<bool>("doSeedPHIHisto");
109  doPHIVsETA = conf_.getParameter<bool>("doSeedPHIVsETAHisto");
110  doTheta = conf_.getParameter<bool>("doSeedThetaHisto");
111  doQ = conf_.getParameter<bool>("doSeedQHisto");
112  doDxy = conf_.getParameter<bool>("doSeedDxyHisto");
113  doDz = conf_.getParameter<bool>("doSeedDzHisto");
114  doNRecHits = conf_.getParameter<bool>("doSeedNRecHitsHisto");
115  doProfPHI = conf_.getParameter<bool>("doSeedNVsPhiProf");
116  doProfETA = conf_.getParameter<bool>("doSeedNVsEtaProf");
117 
118  // if (doAllPlots){doAllSeedPlots=true; doTCPlots=true;}
119 
120  ibooker.setCurrentFolder(MEFolderName);
121 
122  // book the Seed histograms
123  // ---------------------------------------------------------------------------------//
124  // std::cout << "[TrackBuildingAnalyzer::beginRun] MEFolderName: " << MEFolderName << std::endl;
125  ibooker.setCurrentFolder(MEFolderName+"/TrackBuilding");
126 
127  if (doAllSeedPlots || doPT) {
128  histname = "SeedPt_"+seedProducer.label() + "_";
129  SeedPt = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, TrackPtBin, TrackPtMin, TrackPtMax);
130  SeedPt->setAxisTitle("Seed p_{T} (GeV/c)", 1);
131  SeedPt->setAxisTitle("Number of Seeds", 2);
132  }
133 
134  if (doAllSeedPlots || doETA) {
135  histname = "SeedEta_"+seedProducer.label() + "_";
136  SeedEta = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, EtaBin, EtaMin, EtaMax);
137  SeedEta->setAxisTitle("Seed #eta", 1);
138  SeedEta->setAxisTitle("Number of Seeds", 2);
139  }
140 
141  if (doAllSeedPlots || doPHI) {
142  histname = "SeedPhi_"+seedProducer.label() + "_";
143  SeedPhi = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, PhiBin, PhiMin, PhiMax);
144  SeedPhi->setAxisTitle("Seed #phi", 1);
145  SeedPhi->setAxisTitle("Number of Seed", 2);
146  }
147 
148  if (doAllSeedPlots || doPHIVsETA) {
149  histname = "SeedPhiVsEta_"+seedProducer.label() + "_";
150  SeedPhiVsEta = ibooker.book2D(histname+CatagoryName, histname+CatagoryName, EtaBin, EtaMin, EtaMax, PhiBin, PhiMin, PhiMax);
151  SeedPhiVsEta->setAxisTitle("Seed #eta", 1);
152  SeedPhiVsEta->setAxisTitle("Seed #phi", 2);
153  }
154 
155  if (doAllSeedPlots || doTheta){
156  histname = "SeedTheta_"+seedProducer.label() + "_";
157  SeedTheta = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, ThetaBin, ThetaMin, ThetaMax);
158  SeedTheta->setAxisTitle("Seed #theta", 1);
159  SeedTheta->setAxisTitle("Number of Seeds", 2);
160  }
161 
162  if (doAllSeedPlots || doQ){
163  histname = "SeedQ_"+seedProducer.label() + "_";
164  SeedQ = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, TrackQBin, TrackQMin, TrackQMax);
165  SeedQ->setAxisTitle("Seed Charge", 1);
166  SeedQ->setAxisTitle("Number of Seeds",2);
167  }
168 
169  if (doAllSeedPlots || doDxy){
170  histname = "SeedDxy_"+seedProducer.label() + "_";
171  SeedDxy = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, SeedDxyBin, SeedDxyMin, SeedDxyMax);
172  SeedDxy->setAxisTitle("Seed d_{xy} (cm)", 1);
173  SeedDxy->setAxisTitle("Number of Seeds",2);
174  }
175 
176  if (doAllSeedPlots || doDz){
177  histname = "SeedDz_"+seedProducer.label() + "_";
178  SeedDz = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, SeedDzBin, SeedDzMin, SeedDzMax);
179  SeedDz->setAxisTitle("Seed d_{z} (cm)", 1);
180  SeedDz->setAxisTitle("Number of Seeds",2);
181  }
182 
183  if (doAllSeedPlots || doNRecHits){
184  histname = "NumberOfRecHitsPerSeed_"+seedProducer.label() + "_";
185  NumberOfRecHitsPerSeed = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, SeedHitBin, SeedHitMin, SeedHitMax);
186  NumberOfRecHitsPerSeed->setAxisTitle("Number of RecHits per Seed", 1);
187  NumberOfRecHitsPerSeed->setAxisTitle("Number of Seeds",2);
188  }
189 
190  if (doAllSeedPlots || doProfPHI){
191  histname = "NumberOfRecHitsPerSeedVsPhiProfile_"+seedProducer.label() + "_";
192  NumberOfRecHitsPerSeedVsPhiProfile = ibooker.bookProfile(histname+CatagoryName, histname+CatagoryName, PhiBin, PhiMin, PhiMax, SeedHitBin, SeedHitMin, SeedHitMax,"s");
194  NumberOfRecHitsPerSeedVsPhiProfile->setAxisTitle("Number of RecHits of each Seed",2);
195  }
196 
197  if (doAllSeedPlots || doProfETA){
198  histname = "NumberOfRecHitsPerSeedVsEtaProfile_"+seedProducer.label() + "_";
199  NumberOfRecHitsPerSeedVsEtaProfile = ibooker.bookProfile(histname+CatagoryName, histname+CatagoryName, EtaBin, EtaMin, EtaMax, SeedHitBin, SeedHitMin, SeedHitMax,"s");
201  NumberOfRecHitsPerSeedVsEtaProfile->setAxisTitle("Number of RecHits of each Seed",2);
202  }
203 
204  // book the TrackCandidate histograms
205  // ---------------------------------------------------------------------------------//
206 
207  if (doTCPlots){
208 
209  ibooker.setCurrentFolder(MEFolderName+"/TrackBuilding");
210 
211  histname = "TrackCandPt_"+tcProducer.label() + "_";
212  TrackCandPt = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, TrackPtBin, TrackPtMin, TrackPtMax);
213  TrackCandPt->setAxisTitle("Track Candidate p_{T} (GeV/c)", 1);
214  TrackCandPt->setAxisTitle("Number of Track Candidates", 2);
215 
216  histname = "TrackCandEta_"+tcProducer.label() + "_";
217  TrackCandEta = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, EtaBin, EtaMin, EtaMax);
218  TrackCandEta->setAxisTitle("Track Candidate #eta", 1);
219  TrackCandEta->setAxisTitle("Number of Track Candidates", 2);
220 
221  histname = "TrackCandPhi_"+tcProducer.label() + "_";
222  TrackCandPhi = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, PhiBin, PhiMin, PhiMax);
223  TrackCandPhi->setAxisTitle("Track Candidate #phi", 1);
224  TrackCandPhi->setAxisTitle("Number of Track Candidates", 2);
225 
226  if (doTheta) {
227  histname = "TrackCandTheta_"+tcProducer.label() + "_";
228  TrackCandTheta = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, ThetaBin, ThetaMin, ThetaMax);
229  TrackCandTheta->setAxisTitle("Track Candidate #theta", 1);
230  TrackCandTheta->setAxisTitle("Number of Track Candidates", 2);
231  }
232 
233  if (doAllTCPlots) {
234  histname = "TrackCandQ_"+tcProducer.label() + "_";
235  TrackCandQ = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, TrackQBin, TrackQMin, TrackQMax);
236  TrackCandQ->setAxisTitle("Track Candidate Charge", 1);
237  TrackCandQ->setAxisTitle("Number of Track Candidates",2);
238 
239  histname = "TrackCandDxy_"+tcProducer.label() + "_";
240  TrackCandDxy = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, TCDxyBin, TCDxyMin, TCDxyMax);
241  TrackCandDxy->setAxisTitle("Track Candidate d_{xy} (cm)", 1);
242  TrackCandDxy->setAxisTitle("Number of Track Candidates",2);
243 
244  histname = "TrackCandDz_"+tcProducer.label() + "_";
245  TrackCandDz = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, TCDzBin, TCDzMin, TCDzMax);
246  TrackCandDz->setAxisTitle("Track Candidate d_{z} (cm)", 1);
247  TrackCandDz->setAxisTitle("Number of Track Candidates",2);
248 
249  histname = "NumberOfRecHitsPerTrackCand_"+tcProducer.label() + "_";
250  NumberOfRecHitsPerTrackCand = ibooker.book1D(histname+CatagoryName, histname+CatagoryName, TCHitBin, TCHitMin, TCHitMax);
251  NumberOfRecHitsPerTrackCand->setAxisTitle("Number of RecHits per Track Candidate", 1);
252  NumberOfRecHitsPerTrackCand->setAxisTitle("Number of Track Candidates",2);
253 
254  histname = "NumberOfRecHitsPerTrackCandVsPhiProfile_"+tcProducer.label() + "_";
255  NumberOfRecHitsPerTrackCandVsPhiProfile = ibooker.bookProfile(histname+CatagoryName, histname+CatagoryName, PhiBin, PhiMin, PhiMax, TCHitBin, TCHitMin, TCHitMax,"s");
256  NumberOfRecHitsPerTrackCandVsPhiProfile->setAxisTitle("Track Candidate #phi",1);
257  NumberOfRecHitsPerTrackCandVsPhiProfile->setAxisTitle("Number of RecHits of each Track Candidate",2);
258 
259  histname = "NumberOfRecHitsPerTrackCandVsEtaProfile_"+tcProducer.label() + "_";
260  NumberOfRecHitsPerTrackCandVsEtaProfile = ibooker.bookProfile(histname+CatagoryName, histname+CatagoryName, EtaBin, EtaMin, EtaMax, TCHitBin, TCHitMin, TCHitMax,"s");
261  NumberOfRecHitsPerTrackCandVsEtaProfile->setAxisTitle("Track Candidate #eta",1);
262  NumberOfRecHitsPerTrackCandVsEtaProfile->setAxisTitle("Number of RecHits of each Track Candidate",2);
263  }
264 
265  histname = "TrackCandPhiVsEta_"+tcProducer.label() + "_";
266  TrackCandPhiVsEta = ibooker.book2D(histname+CatagoryName, histname+CatagoryName, EtaBin, EtaMin, EtaMax, PhiBin, PhiMin, PhiMax);
267  TrackCandPhiVsEta->setAxisTitle("Track Candidate #eta", 1);
268  TrackCandPhiVsEta->setAxisTitle("Track Candidate #phi", 2);
269  }
270 
271 }
272 
273 // -- Analyse
274 // ---------------------------------------------------------------------------------//
276 (
277  const edm::Event& iEvent,
278  const edm::EventSetup& iSetup,
279  const TrajectorySeed& candidate,
280  const reco::BeamSpot& bs,
281  const edm::ESHandle<MagneticField>& theMF,
283 )
284 {
285  TSCBLBuilderNoMaterial tscblBuilder;
286 
287  //get parameters and errors from the candidate state
288  auto const & theG = ((TkTransientTrackingRecHitBuilder const *)(theTTRHBuilder.product()))->geometry();
289  auto const & candSS = candidate.startingState();
290  TrajectoryStateOnSurface state = trajectoryStateTransform::transientState( candSS, &(theG->idToDet(candSS.detId())->surface()), theMF.product());
291  TrajectoryStateClosestToBeamLine tsAtClosestApproachSeed = tscblBuilder(*state.freeState(),bs);//as in TrackProducerAlgorithm
292  if(!(tsAtClosestApproachSeed.isValid())) {
293  edm::LogVerbatim("TrackBuilding") << "TrajectoryStateClosestToBeamLine not valid";
294  return;
295  }
296  GlobalPoint v0 = tsAtClosestApproachSeed.trackStateAtPCA().position();
297  GlobalVector p = tsAtClosestApproachSeed.trackStateAtPCA().momentum();
298  GlobalPoint v(v0.x()-bs.x0(),v0.y()-bs.y0(),v0.z()-bs.z0());
299 
300  double pt = sqrt(state.globalMomentum().perp2());
301  double eta = state.globalPosition().eta();
302  double phi = state.globalPosition().phi();
303  double theta = state.globalPosition().theta();
304  //double pm = sqrt(state.globalMomentum().mag2());
305  //double pz = state.globalMomentum().z();
306  //double qoverp = tsAtClosestApproachSeed.trackStateAtPCA().charge()/p.mag();
307  //double theta = p.theta();
308  //double lambda = M_PI/2-p.theta();
309  double numberOfHits = candidate.recHits().second-candidate.recHits().first;
310  double dxy = (-v.x()*sin(p.phi())+v.y()*cos(p.phi()));
311  double dz = v.z() - (v.x()*p.x()+v.y()*p.y())/p.perp() * p.z()/p.perp();
312 
313  // fill the ME's
314  if (doAllSeedPlots || doQ)SeedQ->Fill( state.charge() );
315  if (doAllSeedPlots || doPT) SeedPt->Fill( pt );
316  if (doAllSeedPlots || doETA) SeedEta->Fill( eta );
317  if (doAllSeedPlots || doPHI) SeedPhi->Fill( phi );
318  if (doAllSeedPlots || doPHIVsETA) SeedPhiVsEta->Fill( eta, phi);
319  if (doAllSeedPlots || doTheta) SeedTheta->Fill( theta );
320  if (doAllSeedPlots || doDxy) SeedDxy->Fill( dxy );
321  if (doAllSeedPlots || doDz) SeedDz->Fill( dz );
322  if (doAllSeedPlots || doNRecHits) NumberOfRecHitsPerSeed->Fill( numberOfHits );
323  if (doAllSeedPlots || doProfETA) NumberOfRecHitsPerSeedVsEtaProfile->Fill( eta, numberOfHits );
324  if (doAllSeedPlots || doProfPHI) NumberOfRecHitsPerSeedVsPhiProfile->Fill( phi, numberOfHits );
325 
326 }
327 
328 // -- Analyse
329 // ---------------------------------------------------------------------------------//
331 (
332  const edm::Event& iEvent,
333  const edm::EventSetup& iSetup,
334  const TrackCandidate& candidate,
335  const reco::BeamSpot& bs,
336  const edm::ESHandle<MagneticField>& theMF,
338 )
339 {
340  TSCBLBuilderNoMaterial tscblBuilder;
341 
342  //get parameters and errors from the candidate state
343  auto const & theG = ((TkTransientTrackingRecHitBuilder const *)(theTTRHBuilder.product()))->geometry();
344  auto const & candSS = candidate.trajectoryStateOnDet();
345  TrajectoryStateOnSurface state = trajectoryStateTransform::transientState( candSS, &(theG->idToDet(candSS.detId())->surface()), theMF.product());
346  TrajectoryStateClosestToBeamLine tsAtClosestApproachTrackCand = tscblBuilder(*state.freeState(),bs);//as in TrackProducerAlgorithm
347  if(!(tsAtClosestApproachTrackCand.isValid())) {
348  edm::LogVerbatim("TrackBuilding") << "TrajectoryStateClosestToBeamLine not valid";
349  return;
350  }
351  GlobalPoint v0 = tsAtClosestApproachTrackCand.trackStateAtPCA().position();
352  GlobalVector p = tsAtClosestApproachTrackCand.trackStateAtPCA().momentum();
353  GlobalPoint v(v0.x()-bs.x0(),v0.y()-bs.y0(),v0.z()-bs.z0());
354 
355  double pt = sqrt(state.globalMomentum().perp2());
356  double eta = state.globalPosition().eta();
357  double phi = state.globalPosition().phi();
358  double theta = state.globalPosition().theta();
359  //double pm = sqrt(state.globalMomentum().mag2());
360  //double pz = state.globalMomentum().z();
361  //double qoverp = tsAtClosestApproachTrackCand.trackStateAtPCA().charge()/p.mag();
362  //double theta = p.theta();
363  //double lambda = M_PI/2-p.theta();
364  double numberOfHits = candidate.recHits().second-candidate.recHits().first;
365  double dxy = (-v.x()*sin(p.phi())+v.y()*cos(p.phi()));
366 
367  double dz = v.z() - (v.x()*p.x()+v.y()*p.y())/p.perp() * p.z()/p.perp();
368 
369  if (doTCPlots){
370  // fill the ME's
371  if (doAllTCPlots) TrackCandQ->Fill( state.charge() );
372  TrackCandPt->Fill( pt );
373  TrackCandEta->Fill( eta );
374  TrackCandPhi->Fill( phi );
375  TrackCandPhiVsEta->Fill( eta, phi );
376  if (doTheta) TrackCandTheta->Fill( theta );
377  if (doAllTCPlots) TrackCandDxy->Fill( dxy );
378  if (doAllTCPlots) TrackCandDz->Fill( dz );
379  if (doAllTCPlots) NumberOfRecHitsPerTrackCand->Fill( numberOfHits );
380  if (doAllTCPlots) NumberOfRecHitsPerTrackCandVsEtaProfile->Fill( eta, numberOfHits );
381  if (doAllTCPlots) NumberOfRecHitsPerTrackCandVsPhiProfile->Fill( phi, numberOfHits );
382  }
383 }
MonitorElement * NumberOfRecHitsPerTrackCandVsPhiProfile
MonitorElement * TrackCandPt
T getParameter(std::string const &) const
double z0() const
z coordinate
Definition: BeamSpot.h:68
T perp() const
Definition: PV3DBase.h:72
range recHits() const
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:157
MonitorElement * TrackCandDz
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
Geom::Theta< T > theta() const
T y() const
Definition: PV3DBase.h:63
GlobalPoint globalPosition() const
#define NULL
Definition: scimark2.h:8
T perp2() const
Definition: PV3DBase.h:71
MonitorElement * TrackCandPhi
MonitorElement * TrackCandPhiVsEta
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
MonitorElement * NumberOfRecHitsPerTrackCandVsEtaProfile
PTrajectoryStateOnDet const & trajectoryStateOnDet() const
MonitorElement * NumberOfRecHitsPerSeed
virtual void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup, const TrajectorySeed &seed, const reco::BeamSpot &bs, const edm::ESHandle< MagneticField > &theMF, const edm::ESHandle< TransientTrackingRecHitBuilder > &theTTRHBuilder)
int iEvent
Definition: GenABIO.cc:230
MonitorElement * NumberOfRecHitsPerSeedVsPhiProfile
MonitorElement * TrackCandDxy
T sqrt(T t)
Definition: SSEVec.h:48
FreeTrajectoryState const * freeState(bool withErrors=true) const
T z() const
Definition: PV3DBase.h:64
virtual void initHisto(DQMStore::IBooker &ibooker)
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
MonitorElement * TrackCandEta
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
const double EtaMin[kNumberCalorimeter]
TrackBuildingAnalyzer(const edm::ParameterSet &)
PTrajectoryStateOnDet const & startingState() const
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
TrajectoryStateOnSurface transientState(const PTrajectoryStateOnDet &ts, const Surface *surface, const MagneticField *field)
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
T const * product() const
Definition: ESHandle.h:86
range recHits() const
MonitorElement * TrackCandTheta
MonitorElement * SeedPhiVsEta
Geom::Phi< T > phi() const
T eta() const
Definition: PV3DBase.h:76
GlobalVector globalMomentum() const
double y0() const
y coordinate
Definition: BeamSpot.h:66
MonitorElement * NumberOfRecHitsPerSeedVsEtaProfile
T x() const
Definition: PV3DBase.h:62
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
MonitorElement * NumberOfRecHitsPerTrackCand
double x0() const
x coordinate
Definition: BeamSpot.h:64