CMS 3D CMS Logo

LhcTrackAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: LhcTrackAnalyzer
4 // Class: LhcTrackAnalyzer
5 //
13 //
14 
15 // system include files
16 #include <memory>
17 #include <string>
18 #include <vector>
19 
20 // user include files
34 
35 // ROOT includes
36 #include "TFile.h"
37 #include "TTree.h"
38 
39 //
40 // class decleration
41 //
42 
44 public:
45  explicit LhcTrackAnalyzer(const edm::ParameterSet&);
46  ~LhcTrackAnalyzer() override = default;
47  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
48 
49 private:
50  void beginJob() override;
51  void analyze(const edm::Event&, const edm::EventSetup&) override;
52  void endJob() override;
53 
54  // ----------member data ---------------------------
57  bool debug_;
58  std::vector<unsigned int> acceptedBX_;
59 
62 
63  // Output
65  TFile* rootFile_;
66  TTree* rootTree_;
67 
68  // Root-Tuple variables :
69  //=======================
70  void SetVarToZero();
71 
72  static constexpr int nMaxtracks_ = 3000;
73  int nTracks_;
74  int run_;
75  int event_;
76  double pt_[nMaxtracks_];
77  double eta_[nMaxtracks_];
78  double phi_[nMaxtracks_];
79  double chi2_[nMaxtracks_];
83  double dz_[nMaxtracks_];
84  double dxy_[nMaxtracks_];
85  double xPCA_[nMaxtracks_];
86  double yPCA_[nMaxtracks_];
87  double zPCA_[nMaxtracks_];
92  bool goodbx_;
93  bool goodvtx_;
94 };
95 
96 // Constructor
97 
99  : TrackCollectionTag_(iConfig.getParameter<edm::InputTag>("TrackCollectionTag")),
100  PVtxCollectionTag_(iConfig.getParameter<edm::InputTag>("PVtxCollectionTag")),
101  debug_(iConfig.getParameter<bool>("Debug")),
102  acceptedBX_(iConfig.getParameter<std::vector<unsigned int>>("acceptedBX")),
103  filename_(iConfig.getParameter<std::string>("OutputFileName")) {
104  //now do what ever initialization is needed
105  theTrackCollectionToken = consumes<reco::TrackCollection>(TrackCollectionTag_);
106  theVertexCollectionToken = consumes<reco::VertexCollection>(PVtxCollectionTag_);
107 }
108 
109 //
110 // member functions
111 //
112 
113 /*****************************************************************************/
115 /*****************************************************************************/
116 {
117  using namespace edm;
118  using namespace reco;
119  using namespace std;
120 
121  //=======================================================
122  // Initialize Root-tuple variables
123  //=======================================================
124 
125  SetVarToZero();
126 
127  //=======================================================
128  // Retrieve the Track information
129  //=======================================================
130 
131  const auto& vertices = iEvent.get(theVertexCollectionToken);
132  const auto& vtx = vertices.front();
133  if (vtx.isFake()) {
134  goodvtx_ = false;
135  } else {
136  goodvtx_ = true;
137  }
138 
139  int bx = iEvent.bunchCrossing();
140  if (acceptedBX_.empty()) {
141  goodbx_ = true;
142  } else {
143  if (std::find(acceptedBX_.begin(), acceptedBX_.end(), bx) != acceptedBX_.end()) {
144  goodbx_ = true;
145  }
146  }
147 
148  run_ = iEvent.id().run();
149  event_ = iEvent.id().event();
150 
151  const auto& tracks = iEvent.get(theTrackCollectionToken);
152  if (debug_) {
153  edm::LogInfo("LhcTrackAnalyzer") << "LhcTrackAnalyzer::analyze() looping over " << tracks.size() << "tracks."
154  << endl;
155  }
156 
157  for (const auto& track : tracks) {
158  if (nTracks_ >= nMaxtracks_) {
159  edm::LogWarning("LhcTrackAnalyzer")
160  << " LhcTrackAnalyzer::analyze() : Warning - Run " << run_ << " Event " << event_
161  << "\tNumber of tracks: " << tracks.size() << " , greater than " << nMaxtracks_ << std::endl;
162  continue;
163  }
164  pt_[nTracks_] = track.pt();
165  eta_[nTracks_] = track.eta();
166  phi_[nTracks_] = track.phi();
167  chi2_[nTracks_] = track.chi2();
168  chi2ndof_[nTracks_] = track.normalizedChi2();
169  charge_[nTracks_] = track.charge();
170  qoverp_[nTracks_] = track.qoverp();
171  dz_[nTracks_] = track.dz();
172  dxy_[nTracks_] = track.dxy();
173  xPCA_[nTracks_] = track.vertex().x();
174  yPCA_[nTracks_] = track.vertex().y();
175  zPCA_[nTracks_] = track.vertex().z();
176  validhits_[nTracks_][0] = track.numberOfValidHits();
177  validhits_[nTracks_][1] = track.hitPattern().numberOfValidPixelBarrelHits();
178  validhits_[nTracks_][2] = track.hitPattern().numberOfValidPixelEndcapHits();
179  validhits_[nTracks_][3] = track.hitPattern().numberOfValidStripTIBHits();
180  validhits_[nTracks_][4] = track.hitPattern().numberOfValidStripTIDHits();
181  validhits_[nTracks_][5] = track.hitPattern().numberOfValidStripTOBHits();
182  validhits_[nTracks_][6] = track.hitPattern().numberOfValidStripTECHits();
183 
184  int myalgo = -88;
185  if (track.algo() == reco::TrackBase::undefAlgorithm) {
186  myalgo = 0;
187  } else if (track.algo() == reco::TrackBase::ctf) {
188  myalgo = 1;
189  } else if (track.algo() == reco::TrackBase::duplicateMerge) {
190  myalgo = 2;
191  } else if (track.algo() == reco::TrackBase::cosmics) {
192  myalgo = 3;
193  } else if (track.algo() == reco::TrackBase::initialStep) {
194  myalgo = 4;
195  } else if (track.algo() == reco::TrackBase::lowPtTripletStep) {
196  myalgo = 5;
197  } else if (track.algo() == reco::TrackBase::pixelPairStep) {
198  myalgo = 6;
199  } else if (track.algo() == reco::TrackBase::detachedTripletStep) {
200  myalgo = 7;
201  } else if (track.algo() == reco::TrackBase::mixedTripletStep) {
202  myalgo = 8;
203  } else if (track.algo() == reco::TrackBase::pixelLessStep) {
204  myalgo = 9;
205  } else if (track.algo() == reco::TrackBase::tobTecStep) {
206  myalgo = 10;
207  } else if (track.algo() == reco::TrackBase::jetCoreRegionalStep) {
208  myalgo = 11;
209  } else if (track.algo() == reco::TrackBase::muonSeededStepInOut) {
210  myalgo = 13;
211  } else if (track.algo() == reco::TrackBase::muonSeededStepOutIn) {
212  myalgo = 14;
213  } else if (track.algo() == reco::TrackBase::highPtTripletStep) {
214  myalgo = 22;
215  } else if (track.algo() == reco::TrackBase::lowPtQuadStep) {
216  myalgo = 23;
217  } else if (track.algo() == reco::TrackBase::detachedQuadStep) {
218  myalgo = 24;
219  } else {
220  myalgo = 25;
221  edm::LogWarning("LhcTrackAnalyzer")
222  << "LhcTrackAnalyzer does not support all types of tracks, encountered one from algo "
223  << reco::TrackBase::algoName(track.algo());
224  }
225  trkAlgo_[nTracks_] = myalgo;
226 
227  int myquality = -99;
229  myquality = -1;
230  if (track.quality(reco::TrackBase::loose))
231  myquality = 0;
232  if (track.quality(reco::TrackBase::tight))
233  myquality = 1;
234  if (track.quality(reco::TrackBase::highPurity))
235  myquality = 2;
236  trkQuality_[nTracks_] = myquality;
237 
238  if (track.quality(reco::TrackBase::highPurity))
239  isHighPurity_[nTracks_] = 1;
240  else
241  isHighPurity_[nTracks_] = 0;
242  nTracks_++;
243 
244  } //end loop on tracks
245 
246  for (int d = 0; d < nTracks_; ++d) {
247  if (abs(trkQuality_[d]) > 5)
248  edm::LogInfo("LhcTrackAnalyzer") << "MYQUALITY!!! " << trkQuality_[d] << " at track # " << d << "/" << nTracks_
249  << endl;
250  }
251 
252  rootTree_->Fill();
253 }
254 
255 /*****************************************************************************/
257 /*****************************************************************************/
258 {
259  edm::LogInfo("beginJob") << "Begin Job" << std::endl;
260  // Define TTree for output
261  rootFile_ = new TFile(filename_.c_str(), "recreate");
262  rootTree_ = new TTree("tree", "Lhc Track tree");
263 
264  // Track Paramters
265  rootTree_->Branch("run", &run_, "run/I");
266  rootTree_->Branch("event", &event_, "event/I");
267  rootTree_->Branch("goodbx", &goodbx_, "goodbx/O");
268  rootTree_->Branch("goodvtx", &goodvtx_, "goodvtx/O");
269  rootTree_->Branch("nTracks", &nTracks_, "nTracks/I");
270  rootTree_->Branch("pt", &pt_, "pt[nTracks]/D");
271  rootTree_->Branch("eta", &eta_, "eta[nTracks]/D");
272  rootTree_->Branch("phi", &phi_, "phi[nTracks]/D");
273  rootTree_->Branch("chi2", &chi2_, "chi2[nTracks]/D");
274  rootTree_->Branch("chi2ndof", &chi2ndof_, "chi2ndof[nTracks]/D");
275  rootTree_->Branch("charge", &charge_, "charge[nTracks]/I");
276  rootTree_->Branch("qoverp", &qoverp_, "qoverp[nTracks]/D");
277  rootTree_->Branch("dz", &dz_, "dz[nTracks]/D");
278  rootTree_->Branch("dxy", &dxy_, "dxy[nTracks]/D");
279  rootTree_->Branch("xPCA", &xPCA_, "xPCA[nTracks]/D");
280  rootTree_->Branch("yPCA", &yPCA_, "yPCA[nTracks]/D");
281  rootTree_->Branch("zPCA", &zPCA_, "zPCA[nTracks]/D");
282  rootTree_->Branch("isHighPurity", &isHighPurity_, "isHighPurity[nTracks]/I");
283  rootTree_->Branch("trkQuality", &trkQuality_, "trkQuality[nTracks]/I");
284  rootTree_->Branch("trkAlgo", &trkAlgo_, "trkAlgo[nTracks]/I");
285  rootTree_->Branch("nValidHits", &validhits_, "nValidHits[nTracks][7]/I");
286 }
287 
288 /*****************************************************************************/
290 /*****************************************************************************/
291 {
292  if (rootFile_) {
293  rootFile_->Write();
294  rootFile_->Close();
295  }
296 }
297 
298 /*****************************************************************************/
300 /*****************************************************************************/
301 {
303  desc.setComment("Ntuplizer for LHC tracks");
304  desc.add<edm::InputTag>("TrackCollectionTag", edm::InputTag("ALCARECOTkAlMinBias"));
305  desc.add<edm::InputTag>("PVtxCollectionTag", edm::InputTag("offlinePrimaryVertices"));
306  desc.add<bool>("Debug", false);
307  desc.add<std::vector<unsigned int>>("acceptedBX", {});
308  desc.add<std::string>("OutputFileName", "LhcTrackAnalyzer_Output_default.root");
309  descriptions.addWithDefaultLabel(desc);
310 }
311 
312 /*****************************************************************************/
314 /*****************************************************************************/
315 {
316  run_ = -1;
317  event_ = -99;
318  nTracks_ = 0;
319  for (int i = 0; i < nMaxtracks_; ++i) {
320  pt_[i] = 0;
321  eta_[i] = 0;
322  phi_[i] = 0;
323  chi2_[i] = 0;
324  chi2ndof_[i] = 0;
325  charge_[i] = 0;
326  qoverp_[i] = 0;
327  dz_[i] = 0;
328  dxy_[i] = 0;
329  xPCA_[i] = 0;
330  yPCA_[i] = 0;
331  zPCA_[i] = 0;
332  trkQuality_[i] = 0;
333  trkAlgo_[i] = -1;
334  isHighPurity_[i] = -3;
335  for (int j = 0; j < 7; j++) {
336  validhits_[nTracks_][j] = -1 * j;
337  }
338  }
339 }
340 
341 //define this as a plug-in
edm::EDGetTokenT< reco::TrackCollection > theTrackCollectionToken
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
int validhits_[nMaxtracks_][7]
int isHighPurity_[nMaxtracks_]
double dz_[nMaxtracks_]
double zPCA_[nMaxtracks_]
edm::InputTag TrackCollectionTag_
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
int trkQuality_[nMaxtracks_]
double qoverp_[nMaxtracks_]
std::vector< unsigned int > acceptedBX_
edm::EDGetTokenT< reco::VertexCollection > theVertexCollectionToken
~LhcTrackAnalyzer() override=default
double eta_[nMaxtracks_]
void endJob() override
int iEvent
Definition: GenABIO.cc:224
std::string filename_
double phi_[nMaxtracks_]
double chi2ndof_[nMaxtracks_]
double yPCA_[nMaxtracks_]
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static constexpr int nMaxtracks_
double dxy_[nMaxtracks_]
d
Definition: ztail.py:151
std::string algoName() const
Definition: TrackBase.h:550
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Log< level::Info, false > LogInfo
int trkAlgo_[nMaxtracks_]
edm::InputTag PVtxCollectionTag_
double xPCA_[nMaxtracks_]
fixed size matrix
HLT enums.
double chi2_[nMaxtracks_]
int charge_[nMaxtracks_]
Log< level::Warning, false > LogWarning
LhcTrackAnalyzer(const edm::ParameterSet &)
double pt_[nMaxtracks_]
void beginJob() override
void analyze(const edm::Event &, const edm::EventSetup &) override