CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Private Attributes
PatTrackAnalyzer Class Reference
Inheritance diagram for PatTrackAnalyzer:
edm::EDAnalyzer

Classes

struct  Plots
 

Public Member Functions

virtual void analyze (const edm::Event &event, const edm::EventSetup &es)
 
virtual void beginJob ()
 
 PatTrackAnalyzer (const edm::ParameterSet &params)
 constructor and destructor More...
 
 ~PatTrackAnalyzer ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Attributes

edm::InputTag beamSpot_
 
std::vector< Plotsplots_
 
std::vector< std::string > qualities_
 
edm::InputTag src_
 

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)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Definition at line 25 of file PatTrackAnalyzer.cc.

Constructor & Destructor Documentation

PatTrackAnalyzer::PatTrackAnalyzer ( const edm::ParameterSet params)

constructor and destructor

Definition at line 60 of file PatTrackAnalyzer.cc.

60  :
61  src_(params.getParameter<edm::InputTag>("src")),
62  beamSpot_(params.getParameter<edm::InputTag>("beamSpot")),
63  qualities_(params.getParameter< std::vector<std::string> >("qualities"))
64 {
65 }
T getParameter(std::string const &) const
std::vector< std::string > qualities_
edm::InputTag beamSpot_
edm::InputTag src_
PatTrackAnalyzer::~PatTrackAnalyzer ( )

Definition at line 67 of file PatTrackAnalyzer.cc.

68 {
69 }

Member Function Documentation

void PatTrackAnalyzer::analyze ( const edm::Event event,
const edm::EventSetup es 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 138 of file PatTrackAnalyzer.cc.

References abs, beamSpot_, PatTrackAnalyzer::Plots::d0, PatTrackAnalyzer::Plots::d0Err, reco::TrackBase::dxy(), reco::TrackBase::dxyError(), PatTrackAnalyzer::Plots::eta, reco::TrackBase::eta(), reco::TrackBase::hitPattern(), i, PatTrackAnalyzer::Plots::invPt, PatTrackAnalyzer::Plots::invPtErr, edm::Ref< C, T, F >::isNonnull(), edm::HandleBase::isValid(), PatTrackAnalyzer::Plots::nHits, reco::TrackBase::numberOfValidHits(), reco::HitPattern::numberOfValidPixelBarrelHits(), reco::HitPattern::numberOfValidPixelEndcapHits(), reco::HitPattern::numberOfValidStripTECHits(), reco::HitPattern::numberOfValidStripTIBHits(), reco::HitPattern::numberOfValidStripTIDHits(), reco::HitPattern::numberOfValidStripTOBHits(), PatTrackAnalyzer::Plots::phi, reco::TrackBase::phi(), plots_, PatTrackAnalyzer::Plots::pt, reco::TrackBase::pt(), PatTrackAnalyzer::Plots::ptErr, reco::TrackBase::ptError(), PatTrackAnalyzer::Plots::pxbHitsEta, PatTrackAnalyzer::Plots::pxeHitsEta, reco::TrackBase::qoverp(), reco::TrackBase::qoverpError(), qualities_, reco::TrackBase::quality(), reco::TrackBase::qualityByName(), src_, PatTrackAnalyzer::Plots::tecHitsEta, PatTrackAnalyzer::Plots::tibHitsEta, PatTrackAnalyzer::Plots::tidHitsEta, PatTrackAnalyzer::Plots::tobHitsEta, and testEve_cfg::tracks.

139 {
140  // handles to kinds of data we might want to read
144 
145  // read the beam spot
146  event.getByLabel(beamSpot_, beamSpot);
147 
148  // our internal copy of track points
149  // (we need this in order to able to simultaneously access tracks
150  // directly or embedded in PAT objects, like muons, normally you
151  // would iterate over the handle directly)
152  std::vector<const reco::Track*> tracks;
153 
154  event.getByLabel(src_, tracksHandle);
155  if (tracksHandle.isValid()) {
156  // framework was able to read the collection as a view of
157  // tracks, no copy them to our "tracks" variable
158  for(edm::View<reco::Track>::const_iterator iter = tracksHandle->begin();
159  iter != tracksHandle->end(); ++iter)
160  tracks.push_back(&*iter);
161  } else {
162  // does not exist or is not a track collection
163  // let's assume it is a collection of PAT muons
164  event.getByLabel(src_, muonsHandle);
165 
166  // and copy them over
167  // NOTE: We are using ->globalTrack() here
168  // This means we are using the global fit over both
169  // the inner tracker and the muon stations!
170  // other alternatives are: innerTrack(), outerTrack()
171  for(pat::MuonCollection::const_iterator iter = muonsHandle->begin();
172  iter != muonsHandle->end(); ++iter) {
173  reco::TrackRef track = iter->globalTrack();
174  // the muon might not be a "global" muon
175  if (track.isNonnull())
176  tracks.push_back(&*track);
177  }
178  }
179 
180  // we are done filling the tracks into our "tracks" vector.
181  // now analyze them, once for each track quality category
182 
183  unsigned int nQualities = qualities_.size();
184  for(unsigned int i = 0; i < nQualities; ++i) {
185  // we convert the quality flag from its name as a string
186  // to the enumeration value used by the tracking code
187  // (which is essentially an integer number)
189 
190  // our set of plots
191  Plots &plots = plots_[i];
192 
193  // now loop over the tracks
194  for(std::vector<const reco::Track*>::const_iterator iter = tracks.begin();
195  iter != tracks.end(); ++iter) {
196  // this is our track
197  const reco::Track &track = **iter;
198 
199  // ignore tracks that fail the quality cut
200  if (!track.quality(quality))
201  continue;
202 
203  // and fill all the plots
204  plots.eta->Fill(track.eta());
205  plots.phi->Fill(track.phi());
206 
207  plots.pt->Fill(track.pt());
208  plots.ptErr->Fill(track.ptError());
209 
210  plots.invPt->Fill(track.qoverp());
211  plots.invPtErr->Fill(track.qoverpError());
212 
213  // the transverse IP is taken with respect to
214  // the beam spot instead of (0, 0)
215  // because the beam spot in CMS is not at (0, 0)
216  plots.d0->Fill(track.dxy(beamSpot->position()));
217  plots.d0Err->Fill(track.dxyError());
218 
219  plots.nHits->Fill(track.numberOfValidHits());
220 
221  // the hit pattern contains information about
222  // which modules of the detector have been hit
223  const reco::HitPattern &hits = track.hitPattern();
224 
225  double absEta = std::abs(track.eta());
226  // now fill the number of hits in a layer depending on eta
227  plots.pxbHitsEta->Fill(absEta, hits.numberOfValidPixelBarrelHits());
228  plots.pxeHitsEta->Fill(absEta, hits.numberOfValidPixelEndcapHits());
229  plots.tibHitsEta->Fill(absEta, hits.numberOfValidStripTIBHits());
230  plots.tobHitsEta->Fill(absEta, hits.numberOfValidStripTOBHits());
231  plots.tidHitsEta->Fill(absEta, hits.numberOfValidStripTIDHits());
232  plots.tecHitsEta->Fill(absEta, hits.numberOfValidStripTECHits());
233  }
234  }
235 }
double qoverp() const
q/p
Definition: TrackBase.h:115
int i
Definition: DBlmapReader.cc:9
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
TrackQuality
track quality
Definition: TrackBase.h:95
double dxyError() const
error on dxy
Definition: TrackBase.h:209
std::vector< Plots > plots_
std::vector< std::string > qualities_
#define abs(x)
Definition: mlp_lapack.h:159
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:139
int numberOfValidStripTIDHits() const
Definition: HitPattern.cc:424
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:249
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:141
int numberOfValidPixelEndcapHits() const
Definition: HitPattern.cc:385
double pt() const
track transverse momentum
Definition: TrackBase.h:131
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:194
int numberOfValidStripTECHits() const
Definition: HitPattern.cc:450
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:232
edm::InputTag beamSpot_
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:223
bool isValid() const
Definition: HandleBase.h:76
int numberOfValidStripTIBHits() const
Definition: HitPattern.cc:411
double qoverpError() const
error on signed transverse curvature
Definition: TrackBase.h:192
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:46
tuple tracks
Definition: testEve_cfg.py:39
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:364
int numberOfValidStripTOBHits() const
Definition: HitPattern.cc:437
edm::InputTag src_
int numberOfValidPixelBarrelHits() const
Definition: HitPattern.cc:372
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:121
void PatTrackAnalyzer::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 71 of file PatTrackAnalyzer.cc.

References PatTrackAnalyzer::Plots::d0, PatTrackAnalyzer::Plots::d0Err, PatTrackAnalyzer::Plots::eta, i, PatTrackAnalyzer::Plots::invPt, PatTrackAnalyzer::Plots::invPtErr, M_PI, TFileDirectory::make(), PatTrackAnalyzer::Plots::nHits, PatTrackAnalyzer::Plots::phi, plots_, PatTrackAnalyzer::Plots::pt, PatTrackAnalyzer::Plots::ptErr, PatTrackAnalyzer::Plots::pxbHitsEta, PatTrackAnalyzer::Plots::pxeHitsEta, qualities_, PatTrackAnalyzer::Plots::tecHitsEta, PatTrackAnalyzer::Plots::tibHitsEta, PatTrackAnalyzer::Plots::tidHitsEta, and PatTrackAnalyzer::Plots::tobHitsEta.

72 {
73  // retrieve handle to auxiliary service
74  // used for storing histograms into ROOT file
76 
77  // now book the histograms, for each category
78  unsigned int nQualities = qualities_.size();
79 
80  plots_.resize(nQualities);
81 
82  for(unsigned int i = 0; i < nQualities; ++i) {
83  // the name of the quality flag
84  const char *quality = qualities_[i].c_str();
85 
86  // the set of plots
87  Plots &plots = plots_[i];
88 
89  plots.eta = fs->make<TH1F>(Form("eta_%s", quality),
90  Form("track \\eta (%s)", quality),
91  100, -3, 3);
92  plots.phi = fs->make<TH1F>(Form("phi_%s", quality),
93  Form("track \\phi (%s)", quality),
94  100, -M_PI, +M_PI);
95  plots.pt = fs->make<TH1F>(Form("pt_%s", quality),
96  Form("track p_{T} (%s)", quality),
97  100, 0, 10);
98  plots.ptErr = fs->make<TH1F>(Form("ptErr_%s", quality),
99  Form("track p_{T} error (%s)", quality),
100  100, 0, 1);
101  plots.invPt = fs->make<TH1F>(Form("invPt_%s", quality),
102  Form("track 1/p_{T} (%s)", quality),
103  100, -5, 5);
104  plots.invPtErr = fs->make<TH1F>(Form("invPtErr_%s", quality),
105  Form("track 1/p_{T} error (%s)", quality),
106  100, 0, 0.1);
107  plots.d0 = fs->make<TH1F>(Form("d0_%s", quality),
108  Form("track d0 (%s)", quality),
109  100, 0, 0.1);
110  plots.d0Err = fs->make<TH1F>(Form("d0Err_%s", quality),
111  Form("track d0 error (%s)", quality),
112  100, 0, 0.1);
113  plots.nHits = fs->make<TH1F>(Form("nHits_%s", quality),
114  Form("track number of total hits (%s)", quality),
115  60, 0, 60);
116 
117  plots.pxbHitsEta = fs->make<TProfile>(Form("pxbHitsEta_%s", quality),
118  Form("#hits in Pixel Barrel (%s)", quality),
119  100, 0, 3);
120  plots.pxeHitsEta = fs->make<TProfile>(Form("pxeHitsEta_%s", quality),
121  Form("#hits in Pixel Endcap (%s)", quality),
122  100, 0, 3);
123  plots.tibHitsEta = fs->make<TProfile>(Form("tibHitsEta_%s", quality),
124  Form("#hits in Tracker Inner Barrel (%s)", quality),
125  100, 0, 3);
126  plots.tobHitsEta = fs->make<TProfile>(Form("tobHitsEta_%s", quality),
127  Form("#hits in Tracker Outer Barrel (%s)", quality),
128  100, 0, 3);
129  plots.tidHitsEta = fs->make<TProfile>(Form("tidHitsEta_%s", quality),
130  Form("#hits in Tracker Inner Disk (%s)", quality),
131  100, 0, 3);
132  plots.tecHitsEta = fs->make<TProfile>(Form("tecHitsEta_%s", quality),
133  Form("#hits in Tracker Endcap (%s)", quality),
134  100, 0, 3);
135  }
136 }
int i
Definition: DBlmapReader.cc:9
std::vector< Plots > plots_
std::vector< std::string > qualities_
#define M_PI
Definition: BFit3D.cc:3
T * make() const
make new ROOT object

Member Data Documentation

edm::InputTag PatTrackAnalyzer::beamSpot_
private

Definition at line 38 of file PatTrackAnalyzer.cc.

Referenced by analyze().

std::vector<Plots> PatTrackAnalyzer::plots_
private

Definition at line 56 of file PatTrackAnalyzer.cc.

Referenced by analyze(), and beginJob().

std::vector<std::string> PatTrackAnalyzer::qualities_
private

Definition at line 41 of file PatTrackAnalyzer.cc.

Referenced by analyze(), and beginJob().

edm::InputTag PatTrackAnalyzer::src_
private

Definition at line 37 of file PatTrackAnalyzer.cc.

Referenced by analyze().