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
SegmentTrackAnalyzer Class Reference

#include <SegmentTrackAnalyzer.h>

Inheritance diagram for SegmentTrackAnalyzer:
MuonAnalyzerBase

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &, const reco::Track &recoTrack)
 Get the analysis. More...
 
void beginJob (DQMStore *dbe)
 Inizialize parameters for histo binning. More...
 
 SegmentTrackAnalyzer (const edm::ParameterSet &, MuonServiceProxy *theService)
 Constructor. More...
 
virtual ~SegmentTrackAnalyzer ()
 Destructor. More...
 
- Public Member Functions inherited from MuonAnalyzerBase
void analyze (const edm::Event &, const edm::EventSetup &, reco::Muon &recoMuon)
 Get the analysis of the muon properties. More...
 
void analyze (const edm::Event &, const edm::EventSetup &, reco::Track &recoTrack)
 Get the analysis of the muon track properties. More...
 
 MuonAnalyzerBase (MuonServiceProxy *theServ)
 Constructor. More...
 
MuonServiceProxyservice ()
 
virtual ~MuonAnalyzerBase ()
 Destructor. More...
 

Private Attributes

MonitorElementcscTrackHitPercentualVsEta
 
MonitorElementcscTrackHitPercentualVsPhi
 
MonitorElementcscTrackHitPercentualVsPt
 
MonitorElementdtTrackHitPercentualVsEta
 
MonitorElementdtTrackHitPercentualVsPhi
 
MonitorElementdtTrackHitPercentualVsPt
 
MonitorElementhitsNotUsed
 
MonitorElementhitsNotUsedPercentual
 
MonitorElementhitStaProvenance
 
MonitorElementhitTkrProvenance
 
std::string metname
 
edm::ParameterSet parameters
 
SegmentsTrackAssociatortheSegmentsAssociator
 
MonitorElementtrackHitPercentualVsEta
 
MonitorElementtrackHitPercentualVsPhi
 
MonitorElementtrackHitPercentualVsPt
 
MonitorElementTrackSegm
 

Detailed Description

DQM monitoring source for segments associated to the muon track

Date:
2009/12/22 17:43:41
Revision:
1.8
Author
G. Mila - INFN Torino

Definition at line 29 of file SegmentTrackAnalyzer.h.

Constructor & Destructor Documentation

SegmentTrackAnalyzer::SegmentTrackAnalyzer ( const edm::ParameterSet pSet,
MuonServiceProxy theService 
)

Constructor.

Definition at line 32 of file SegmentTrackAnalyzer.cc.

References edm::ParameterSet::getParameter(), parameters, and theSegmentsAssociator.

33 
34  parameters = pSet;
35 
36  const ParameterSet SegmentsTrackAssociatorParameters = parameters.getParameter<ParameterSet>("SegmentsTrackAssociatorParameters");
37  theSegmentsAssociator = new SegmentsTrackAssociator(SegmentsTrackAssociatorParameters);
38 
39 }
T getParameter(std::string const &) const
MuonServiceProxy * theService
MuonAnalyzerBase(MuonServiceProxy *theServ)
Constructor.
SegmentsTrackAssociator * theSegmentsAssociator
edm::ParameterSet parameters
SegmentTrackAnalyzer::~SegmentTrackAnalyzer ( )
virtual

Destructor.

Definition at line 42 of file SegmentTrackAnalyzer.cc.

42 { }

Member Function Documentation

void SegmentTrackAnalyzer::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup,
const reco::Track recoTrack 
)

Get the analysis.

Definition at line 112 of file SegmentTrackAnalyzer.cc.

References SegmentsTrackAssociator::associate(), MuonSubdetId::CSC, cscTrackHitPercentualVsEta, cscTrackHitPercentualVsPhi, cscTrackHitPercentualVsPt, MuonSubdetId::DT, dtTrackHitPercentualVsEta, dtTrackHitPercentualVsPhi, dtTrackHitPercentualVsPt, reco::TrackBase::eta(), MonitorElement::Fill(), hitsNotUsed, hitsNotUsedPercentual, hitStaProvenance, hitTkrProvenance, LogTrace, metname, DetId::Muon, reco::TrackBase::phi(), DTRecSegment4D::phiSegment(), PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, reco::TrackBase::pt(), reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), MuonSubdetId::RPC, findQualityFiles::size, DTRecSegment2D::specificRecHits(), SiStripDetId::TEC, theSegmentsAssociator, SiStripDetId::TIB, SiStripDetId::TID, SiStripDetId::TOB, DetId::Tracker, trackHitPercentualVsEta, trackHitPercentualVsPhi, trackHitPercentualVsPt, and TrackSegm.

112  {
113 
114  LogTrace(metname)<<"[SegmentTrackAnalyzer] Filling the histos";
115 
117 
118  LogTrace(metname)<<"[SegmentTrackAnalyzer] # of segments associated to the track: "<<(segments).size();
119 
120  // hit counters
121  int hitsFromDt=0;
122  int hitsFromCsc=0;
123  int hitsFromRpc=0;
124  int hitsFromTk=0;
125  int hitsFromTrack=0;
126  int hitsFromSegmDt=0;
127  int hitsFromSegmCsc=0;
128  // segment counters
129  int segmFromDt=0;
130  int segmFromCsc=0;
131 
132  for (MuonTransientTrackingRecHit::MuonRecHitContainer::const_iterator segment=segments.begin();
133  segment!=segments.end(); segment++) {
134 
135  DetId id = (*segment)->geographicalId();
136 
137  // hits from DT segments
138  if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::DT ) {
139  ++segmFromDt;
140  const DTRecSegment4D *seg4D = dynamic_cast<const DTRecSegment4D*>((*segment)->hit());
141  if((*seg4D).hasPhi())
142  hitsFromSegmDt+=(*seg4D).phiSegment()->specificRecHits().size();
143  if((*seg4D).hasZed())
144  hitsFromSegmDt+=(*seg4D).zSegment()->specificRecHits().size();
145 
146  }
147 
148  // hits from CSC segments
149  if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::CSC ) {
150  hitsFromSegmCsc+=(*segment)->recHits().size();
151  segmFromCsc++;
152  }
153 
154  }
155 
156 
157  // hits from track
158  for(trackingRecHit_iterator recHit = recoTrack.recHitsBegin(); recHit != recoTrack.recHitsEnd(); ++recHit){
159 
160  hitsFromTrack++;
161  DetId id = (*recHit)->geographicalId();
162  // hits from DT
163  if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::DT )
164  hitsFromDt++;
165  // hits from CSC
166  if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::CSC )
167  hitsFromCsc++;
168  // hits from RPC
169  if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::RPC )
170  hitsFromRpc++;
171  // hits from Tracker
172  if (id.det() == DetId::Tracker){
173  hitsFromTk++;
174  if(id.subdetId() == PixelSubdetector::PixelBarrel )
176  if(id.subdetId() == PixelSubdetector::PixelEndcap )
178  if(id.subdetId() == SiStripDetId::TIB )
180  if(id.subdetId() == SiStripDetId::TID )
182  if(id.subdetId() == SiStripDetId::TOB )
184  if(id.subdetId() == SiStripDetId::TEC )
186  }
187 
188  }
189 
190  // fill the histos
191  hitsNotUsed->Fill(hitsFromSegmDt+hitsFromSegmCsc+hitsFromRpc+hitsFromTk-hitsFromTrack);
192  hitsNotUsedPercentual->Fill(double(hitsFromSegmDt+hitsFromSegmCsc+hitsFromRpc+hitsFromTk-hitsFromTrack)/hitsFromTrack);
193 
194  if(hitsFromDt!=0 && hitsFromCsc!=0)
195  TrackSegm->Fill(1,segmFromDt+segmFromCsc);
196  if(hitsFromDt!=0 && hitsFromCsc==0)
197  TrackSegm->Fill(2,segmFromDt);
198  if(hitsFromDt==0 && hitsFromCsc!=0)
199  TrackSegm->Fill(3,segmFromCsc);
200 
201  if(hitsFromDt!=0 && hitsFromCsc==0 && hitsFromRpc==0) hitStaProvenance->Fill(1);
202  if(hitsFromCsc!=0 && hitsFromDt==0 && hitsFromRpc==0) hitStaProvenance->Fill(2);
203  if(hitsFromRpc!=0 && hitsFromDt==0 && hitsFromCsc==0) hitStaProvenance->Fill(3);
204  if(hitsFromDt!=0 && hitsFromCsc!=0 && hitsFromRpc==0) hitStaProvenance->Fill(4);
205  if(hitsFromDt!=0 && hitsFromRpc!=0 && hitsFromCsc==0) hitStaProvenance->Fill(5);
206  if(hitsFromCsc!=0 && hitsFromRpc!=0 && hitsFromDt==0) hitStaProvenance->Fill(6);
207  if(hitsFromDt!=0 && hitsFromCsc!=0 && hitsFromRpc!=0) hitStaProvenance->Fill(7);
208 
209  if(hitsFromSegmDt+hitsFromSegmCsc !=0){
210  trackHitPercentualVsEta->Fill(recoTrack.eta(), double(hitsFromDt+hitsFromCsc)/(hitsFromSegmDt+hitsFromSegmCsc));
211  trackHitPercentualVsPhi->Fill(recoTrack.phi(), double(hitsFromDt+hitsFromCsc)/(hitsFromSegmDt+hitsFromSegmCsc));
212  trackHitPercentualVsPt->Fill(recoTrack.pt(), double(hitsFromDt+hitsFromCsc)/(hitsFromSegmDt+hitsFromSegmCsc));
213  }
214 
215  if(hitsFromSegmDt!=0){
216  dtTrackHitPercentualVsEta->Fill(recoTrack.eta(), double(hitsFromDt)/hitsFromSegmDt);
217  dtTrackHitPercentualVsPhi->Fill(recoTrack.phi(), double(hitsFromDt)/hitsFromSegmDt);
218  dtTrackHitPercentualVsPt->Fill(recoTrack.pt(), double(hitsFromDt)/hitsFromSegmDt);
219  }
220 
221  if(hitsFromSegmCsc!=0){
222  cscTrackHitPercentualVsEta->Fill(recoTrack.eta(), double(hitsFromCsc)/hitsFromSegmCsc);
223  cscTrackHitPercentualVsPhi->Fill(recoTrack.phi(), double(hitsFromCsc)/hitsFromSegmCsc);
224  cscTrackHitPercentualVsPt->Fill(recoTrack.pt(), double(hitsFromCsc)/hitsFromSegmCsc);
225  }
226 
227 }
MonitorElement * trackHitPercentualVsPt
const DTChamberRecSegment2D * phiSegment() const
The superPhi segment: 0 if no phi projection available.
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:139
SegmentsTrackAssociator * theSegmentsAssociator
MonitorElement * cscTrackHitPercentualVsPhi
void Fill(long long x)
static const int CSC
Definition: MuonSubdetId.h:15
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:141
double pt() const
track transverse momentum
Definition: TrackBase.h:131
MonitorElement * hitTkrProvenance
MuonTransientTrackingRecHit::MuonRecHitContainer associate(const edm::Event &, const edm::EventSetup &, const reco::Track &)
Get the analysis.
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:63
MonitorElement * trackHitPercentualVsPhi
#define LogTrace(id)
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
Definition: DetId.h:20
MonitorElement * TrackSegm
MonitorElement * cscTrackHitPercentualVsPt
static const int RPC
Definition: MuonSubdetId.h:16
MonitorElement * dtTrackHitPercentualVsPhi
MonitorElement * hitsNotUsedPercentual
static const int DT
Definition: MuonSubdetId.h:14
MonitorElement * cscTrackHitPercentualVsEta
MonitorElement * dtTrackHitPercentualVsEta
MonitorElement * hitsNotUsed
MonitorElement * dtTrackHitPercentualVsPt
MonitorElement * trackHitPercentualVsEta
std::vector< MuonRecHitPointer > MuonRecHitContainer
tuple size
Write out results.
MonitorElement * hitStaProvenance
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:65
void SegmentTrackAnalyzer::beginJob ( DQMStore dbe)
virtual

Inizialize parameters for histo binning.

Implements MuonAnalyzerBase.

Definition at line 45 of file SegmentTrackAnalyzer.cc.

References DQMStore::book1D(), DQMStore::book2D(), cscTrackHitPercentualVsEta, cscTrackHitPercentualVsPhi, cscTrackHitPercentualVsPt, dtTrackHitPercentualVsEta, dtTrackHitPercentualVsPhi, dtTrackHitPercentualVsPt, jptDQMConfig_cff::etaMax, jptDQMConfig_cff::etaMin, edm::ParameterSet::getParameter(), hitsNotUsed, hitsNotUsedPercentual, hitStaProvenance, hitTkrProvenance, instance, diffTwoXMLs::label, LogTrace, metname, parameters, jptDQMConfig_cff::phiMax, jptDQMConfig_cff::phiMin, jptDQMConfig_cff::ptMax, PtMinSelector_cfg::ptMin, MonitorElement::setAxisTitle(), MonitorElement::setBinLabel(), DQMStore::setCurrentFolder(), trackHitPercentualVsEta, trackHitPercentualVsPhi, trackHitPercentualVsPt, and TrackSegm.

45  {
46 
47 
48  metname = "segmTrackAnalyzer";
49  string trackCollection = parameters.getParameter<edm::InputTag>("MuTrackCollection").label() + parameters.getParameter<edm::InputTag>("MuTrackCollection").instance();
50  LogTrace(metname)<<"[SegmentTrackAnalyzer] Parameters initialization";
51  dbe->setCurrentFolder("Muons/SegmentTrackAnalyzer");
52 
53  // histograms initalization
54  hitsNotUsed = dbe->book1D("HitsNotUsedForGlobalTracking_"+trackCollection, "recHits not used for GLB ["+trackCollection+"]", 50, -0.5, 49.5);
55  hitsNotUsedPercentual = dbe->book1D("HitsNotUsedForGlobalTrackingDvHitUsed_"+trackCollection, "(recHits_{notUsedForGLB}) / (recHits_{GLB}) ["+trackCollection+"]", 100, 0, 1.);
56 
57  TrackSegm = dbe->book2D("trackSegments_"+trackCollection, "Number of segments associated to the track ["+trackCollection+"]", 3, 0.5, 3.5, 8, 0, 8);
58  TrackSegm->setBinLabel(1,"DT+CSC",1);
59  TrackSegm->setBinLabel(2,"DT",1);
60  TrackSegm->setBinLabel(3,"CSC",1);
61 
62  hitStaProvenance = dbe->book1D("trackHitStaProvenance_"+trackCollection, "Number of recHits_{STAinTrack} ["+trackCollection+"]", 7, 0.5, 7.5);
64  hitStaProvenance->setBinLabel(2,"CSC");
65  hitStaProvenance->setBinLabel(3,"RPC");
66  hitStaProvenance->setBinLabel(4,"DT+CSC");
67  hitStaProvenance->setBinLabel(5,"DT+RPC");
68  hitStaProvenance->setBinLabel(6,"CSC+RPC");
69  hitStaProvenance->setBinLabel(7,"DT+CSC+RPC");
70 
71 
72  if(trackCollection!="standAloneMuons"){
73  hitTkrProvenance = dbe->book1D("trackHitTkrProvenance_"+trackCollection, "Number of recHits_{TKinTrack} ["+trackCollection+"]", 6, 0.5, 6.5);
74  hitTkrProvenance->setBinLabel(1,"PixBarrel");
75  hitTkrProvenance->setBinLabel(2,"PixEndCap");
76  hitTkrProvenance->setBinLabel(3,"TIB");
77  hitTkrProvenance->setBinLabel(4,"TID");
78  hitTkrProvenance->setBinLabel(5,"TOB");
79  hitTkrProvenance->setBinLabel(6,"TEC");
80  }
81 
82  int etaBin = parameters.getParameter<int>("etaBin");
83  double etaMin = parameters.getParameter<double>("etaMin");
84  double etaMax = parameters.getParameter<double>("etaMax");
85  trackHitPercentualVsEta = dbe->book2D("trackHitDivtrackSegmHitVsEta_"+trackCollection, "(recHits_{Track} / recHits_{associatedSegm}) vs #eta [" +trackCollection+"]", etaBin, etaMin, etaMax, 20, 0, 1);
86  dtTrackHitPercentualVsEta = dbe->book2D("dtTrackHitDivtrackSegmHitVsEta_"+trackCollection, "(recHits_{DTinTrack} / recHits_{associatedSegm}) vs #eta [" +trackCollection+"]", etaBin, etaMin, etaMax, 20, 0, 1);
87  cscTrackHitPercentualVsEta = dbe->book2D("cscTrackHitDivtrackSegmHitVsEta_"+trackCollection, "(recHits_{CSCinTrack} / recHits_{associatedSegm}) vs #eta [" +trackCollection+"]", etaBin, etaMin, etaMax, 20, 0, 1);
88 
89  int phiBin = parameters.getParameter<int>("phiBin");
90  double phiMin = parameters.getParameter<double>("phiMin");
91  double phiMax = parameters.getParameter<double>("phiMax");
92  trackHitPercentualVsPhi = dbe->book2D("trackHitDivtrackSegmHitVsPhi_"+trackCollection, "(recHits_{Track} / recHits_{associatedSegm}) vs #phi [" +trackCollection+"]", phiBin, phiMin, phiMax, 20, 0, 1);
93  trackHitPercentualVsPhi->setAxisTitle("rad",2);
94  dtTrackHitPercentualVsPhi = dbe->book2D("dtTrackHitDivtrackSegmHitVsPhi_"+trackCollection, "(recHits_{DTinTrack} / recHits_{associatedSegm}) vs #phi [" +trackCollection+"]", phiBin, phiMin, phiMax, 20, 0, 1);
96  cscTrackHitPercentualVsPhi = dbe->book2D("cscTrackHitDivtrackSegmHitVsPhi_"+trackCollection, "(recHits_{CSCinTrack} / recHits_{associatedSegm}) vs #phi [" +trackCollection+"]", phiBin, phiMin, phiMax, 20, 0, 1);
98 
99  int ptBin = parameters.getParameter<int>("ptBin");
100  double ptMin = parameters.getParameter<double>("ptMin");
101  double ptMax = parameters.getParameter<double>("ptMax");
102  trackHitPercentualVsPt = dbe->book2D("trackHitDivtrackSegmHitVsPt_"+trackCollection, "(recHits_{Track} / recHits_{associatedSegm}) vs 1/p_{t} [" +trackCollection+"]", ptBin, ptMin, ptMax, 20, 0, 1);
103  trackHitPercentualVsPt->setAxisTitle("GeV",2);
104  dtTrackHitPercentualVsPt = dbe->book2D("dtTrackHitDivtrackSegmHitVsPt_"+trackCollection, "(recHits_{DTinTrack} / recHits_{associatedSegm}) vs 1/p_{t} [" +trackCollection+"]", ptBin, ptMin, ptMax, 20, 0, 1);
106  cscTrackHitPercentualVsPt = dbe->book2D("cscTrackHitDivtrackSegmHitVsPt_"+trackCollection, "(recHits_{CSCinTrack} / recHits_{associatedSegm}) vs 1/p_{t} [" +trackCollection+"]", ptBin, ptMin, ptMax, 20, 0, 1);
108 
109 }
T getParameter(std::string const &) const
MonitorElement * trackHitPercentualVsPt
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:717
static PFTauRenderPlugin instance
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
MonitorElement * cscTrackHitPercentualVsPhi
MonitorElement * hitTkrProvenance
MonitorElement * trackHitPercentualVsPhi
#define LogTrace(id)
MonitorElement * TrackSegm
MonitorElement * cscTrackHitPercentualVsPt
edm::ParameterSet parameters
MonitorElement * dtTrackHitPercentualVsPhi
MonitorElement * hitsNotUsedPercentual
MonitorElement * cscTrackHitPercentualVsEta
MonitorElement * dtTrackHitPercentualVsEta
MonitorElement * hitsNotUsed
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
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
MonitorElement * dtTrackHitPercentualVsPt
MonitorElement * trackHitPercentualVsEta
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
MonitorElement * hitStaProvenance

Member Data Documentation

MonitorElement* SegmentTrackAnalyzer::cscTrackHitPercentualVsEta
private

Definition at line 66 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().

MonitorElement* SegmentTrackAnalyzer::cscTrackHitPercentualVsPhi
private

Definition at line 67 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().

MonitorElement* SegmentTrackAnalyzer::cscTrackHitPercentualVsPt
private

Definition at line 68 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().

MonitorElement* SegmentTrackAnalyzer::dtTrackHitPercentualVsEta
private

Definition at line 63 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().

MonitorElement* SegmentTrackAnalyzer::dtTrackHitPercentualVsPhi
private

Definition at line 64 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().

MonitorElement* SegmentTrackAnalyzer::dtTrackHitPercentualVsPt
private

Definition at line 65 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().

MonitorElement* SegmentTrackAnalyzer::hitsNotUsed
private

Definition at line 55 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().

MonitorElement* SegmentTrackAnalyzer::hitsNotUsedPercentual
private

Definition at line 56 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().

MonitorElement* SegmentTrackAnalyzer::hitStaProvenance
private

Definition at line 58 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().

MonitorElement* SegmentTrackAnalyzer::hitTkrProvenance
private

Definition at line 59 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().

std::string SegmentTrackAnalyzer::metname
private

Definition at line 50 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().

edm::ParameterSet SegmentTrackAnalyzer::parameters
private
SegmentsTrackAssociator* SegmentTrackAnalyzer::theSegmentsAssociator
private

Definition at line 52 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and SegmentTrackAnalyzer().

MonitorElement* SegmentTrackAnalyzer::trackHitPercentualVsEta
private

Definition at line 60 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().

MonitorElement* SegmentTrackAnalyzer::trackHitPercentualVsPhi
private

Definition at line 61 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().

MonitorElement* SegmentTrackAnalyzer::trackHitPercentualVsPt
private

Definition at line 62 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().

MonitorElement* SegmentTrackAnalyzer::TrackSegm
private

Definition at line 57 of file SegmentTrackAnalyzer.h.

Referenced by analyze(), and beginJob().