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