CMS 3D CMS Logo

XtalDedxAnalysis.cc
Go to the documentation of this file.
1 
2 // system include files
3 #include <memory>
4 #include <string>
5 #include <vector>
6 
7 // user include files
10 
15 
19 
21 
24 
31 
32 #include <TH1F.h>
33 #include <TH1I.h>
34 
35 #define EDM_ML_DEBUG
36 
38 public:
39  explicit XtalDedxAnalysis(const edm::ParameterSet &);
40  ~XtalDedxAnalysis() override {}
41 
42 protected:
43  void beginJob() override {}
44  void analyze(const edm::Event &, const edm::EventSetup &) override;
45  void endJob() override {}
46 
47  void analyzeHits(std::vector<PCaloHit> &,
50 
51 private:
54 
58 
59  TH1F *meNHit_[4], *meE1T0_[4], *meE9T0_[4], *meE1T1_[4], *meE9T1_[4];
60  TH1I *mType_;
61 };
62 
64  caloHitSource_ = ps.getParameter<edm::InputTag>("caloHitSource");
65  simTkLabel_ = ps.getUntrackedParameter<std::string>("moduleLabelTk", "g4SimHits");
66  double energyMax = ps.getParameter<double>("EnergyMax");
67 #ifdef EDM_ML_DEBUG
68  edm::LogInfo("CherenkovAnalysis") << "XtalDedxAnalysis::Source " << caloHitSource_ << " Track Label " << simTkLabel_
69  << " Energy Max " << energyMax;
70 #endif
71  // register for data access
72  tok_calo_ = consumes<edm::PCaloHitContainer>(caloHitSource_);
73  tok_tk_ = consumes<edm::SimTrackContainer>(edm::InputTag(simTkLabel_));
74  tok_vtx_ = consumes<edm::SimVertexContainer>(edm::InputTag(simTkLabel_));
75 
76  // Book histograms
78 
79  if (!tfile.isAvailable())
80  throw cms::Exception("BadConfig") << "TFileService unavailable: "
81  << "please add it to config file";
82  // Histograms for Hits
83  std::string types[4] = {"total", "by dE/dx", "by delta-ray", "by bremms"};
84  char name[20], title[80];
85  for (int i = 0; i < 4; i++) {
86  sprintf(name, "Hits%d", i);
87  sprintf(title, "Number of hits (%s)", types[i].c_str());
88  meNHit_[i] = tfile->make<TH1F>(name, title, 5000, 0., 5000.);
89  meNHit_[i]->GetXaxis()->SetTitle(title);
90  meNHit_[i]->GetYaxis()->SetTitle("Events");
91  sprintf(name, "E1T0%d", i);
92  sprintf(title, "E1 (Loss %s) in GeV", types[i].c_str());
93  meE1T0_[i] = tfile->make<TH1F>(name, title, 5000, 0, energyMax);
94  meE1T0_[i]->GetXaxis()->SetTitle(title);
95  meE1T0_[i]->GetYaxis()->SetTitle("Events");
96  sprintf(name, "E9T0%d", i);
97  sprintf(title, "E9 (Loss %s) in GeV", types[i].c_str());
98  meE9T0_[i] = tfile->make<TH1F>(name, title, 5000, 0, energyMax);
99  meE9T0_[i]->GetXaxis()->SetTitle(title);
100  meE9T0_[i]->GetYaxis()->SetTitle("Events");
101  sprintf(name, "E1T1%d", i);
102  sprintf(title, "E1 (Loss %s with t < 400 ns) in GeV", types[i].c_str());
103  meE1T1_[i] = tfile->make<TH1F>(name, title, 5000, 0, energyMax);
104  meE1T1_[i]->GetXaxis()->SetTitle(title);
105  meE1T1_[i]->GetYaxis()->SetTitle("Events");
106  sprintf(name, "E9T1%d", i);
107  sprintf(title, "E9 (Loss %s with t < 400 ns) in GeV", types[i].c_str());
108  meE9T1_[i] = tfile->make<TH1F>(name, title, 5000, 0, energyMax);
109  meE9T1_[i]->GetXaxis()->SetTitle(title);
110  meE9T1_[i]->GetYaxis()->SetTitle("Events");
111  }
112  sprintf(name, "PDGType");
113  sprintf(title, "PDG ID of first level secondary");
114  mType_ = tfile->make<TH1I>(name, title, 5000, -2500, 2500);
115  mType_->GetXaxis()->SetTitle(title);
116  mType_->GetYaxis()->SetTitle("Tracks");
117 }
118 
120 #ifdef EDM_ML_DEBUG
121  edm::LogVerbatim("CherenkovAnalysis") << "XtalDedxAnalysis::Run = " << e.id().run() << " Event = " << e.id().event();
122 #endif
123  std::vector<PCaloHit> caloHits;
125  e.getByToken(tok_calo_, pCaloHits);
126 
127  std::vector<SimTrack> theSimTracks;
129  e.getByToken(tok_tk_, simTk);
130 
131  std::vector<SimVertex> theSimVertex;
133  e.getByToken(tok_vtx_, simVtx);
134 
135  if (pCaloHits.isValid()) {
136  caloHits.insert(caloHits.end(), pCaloHits->begin(), pCaloHits->end());
137 #ifdef EDM_ML_DEBUG
138  edm::LogVerbatim("CherenkovAnalysis") << "XtalDedxAnalysis: Hit buffer " << caloHits.size();
139 #endif
140  analyzeHits(caloHits, simTk, simVtx);
141  }
142 }
143 
144 void XtalDedxAnalysis::analyzeHits(std::vector<PCaloHit> &hits,
147  edm::SimTrackContainer::const_iterator simTrkItr;
148  int nHit = hits.size();
149  double e10[4], e90[4], e11[4], e91[4], hit[4];
150  for (int i = 0; i < 4; i++)
151  e10[i] = e90[i] = e11[i] = e91[i] = hit[i] = 0;
152  for (int i = 0; i < nHit; i++) {
153  double energy = hits[i].energy();
154  double time = hits[i].time();
155  unsigned int id_ = hits[i].id();
156  int trackID = hits[i].geantTrackId();
157  int type = 1;
158  for (simTrkItr = SimTk->begin(); simTrkItr != SimTk->end(); simTrkItr++) {
159  if (trackID == (int)(simTrkItr->trackId())) {
160  int thePID = simTrkItr->type();
161  if (thePID == 11)
162  type = 2;
163  else if (thePID != -13 && thePID != 13)
164  type = 3;
165  break;
166  }
167  }
168  hit[0]++;
169  hit[type]++;
170  e90[0] += energy;
171  e90[type] += energy;
172  if (time < 400) {
173  e91[0] += energy;
174  e91[type] += energy;
175  }
176  if (id_ == 22) {
177  e10[0] += energy;
178  e10[type] += energy;
179  if (time < 400) {
180  e11[0] += energy;
181  e11[type] += energy;
182  }
183  }
184 #ifdef EDM_ML_DEBUG
185  edm::LogVerbatim("CherenkovAnalysis") << "XtalDedxAnalysis:Hit[" << i << "] ID " << id_ << " E " << energy
186  << " time " << time << " track " << trackID << " type " << type;
187 #endif
188  }
189  for (int i = 0; i < 4; i++) {
190 #ifdef EDM_ML_DEBUG
191  edm::LogVerbatim("CherenkovAnalysis") << "XtalDedxAnalysis:Type(" << i << ") Hit " << hit[i] << " E10 " << e10[i]
192  << " E11 " << e11[i] << " E90 " << e90[i] << " E91 " << e91[i];
193 #endif
194  meNHit_[i]->Fill(hit[i]);
195  meE1T0_[i]->Fill(e10[i]);
196  meE9T0_[i]->Fill(e90[i]);
197  meE1T1_[i]->Fill(e11[i]);
198  meE9T1_[i]->Fill(e91[i]);
199  }
200 
201  // Type of the secondary (coming directly from a generator level track)
202  int nvtx = 0, ntrk = 0, k1 = 0;
203  edm::SimVertexContainer::const_iterator simVtxItr;
204  for (simTrkItr = SimTk->begin(); simTrkItr != SimTk->end(); simTrkItr++)
205  ntrk++;
206  for (simVtxItr = SimVtx->begin(); simVtxItr != SimVtx->end(); simVtxItr++)
207  nvtx++;
208 #ifdef EDM_ML_DEBUG
209  edm::LogVerbatim("CherenkovAnalysis") << "XtalDedxAnalysis: " << ntrk << " tracks and " << nvtx << " vertices";
210 #endif
211  for (simTrkItr = SimTk->begin(); simTrkItr != SimTk->end(); simTrkItr++, ++k1) {
212 #ifdef EDM_ML_DEBUG
213  edm::LogVerbatim("CherenkovAnalysis") << "Track " << k1 << " PDGId " << simTrkItr->type() << " Vertex ID "
214  << simTrkItr->vertIndex() << " Generator " << simTrkItr->noGenpart();
215 #endif
216  if (simTrkItr->noGenpart()) { // This is a secondary
217  int vertIndex = simTrkItr->vertIndex(); // Vertex index of origin
218  if (vertIndex >= 0 && vertIndex < nvtx) {
219  simVtxItr = SimVtx->begin();
220  for (int iv = 0; iv < vertIndex; iv++)
221  simVtxItr++;
222  int parent = simVtxItr->parentIndex(), k2 = 0;
223  for (edm::SimTrackContainer::const_iterator trkItr = SimTk->begin(); trkItr != SimTk->end(); trkItr++, ++k2) {
224 #ifdef EDM_ML_DEBUG
225  edm::LogVerbatim("CherenkovAnalysis") << "XtalDedxAnalysis::Track " << k2 << " ID " << trkItr->trackId()
226  << " (" << parent << ") Generator " << trkItr->noGenpart();
227 #endif
228  if ((int)trkItr->trackId() == parent) { // Parent track
229  if (!trkItr->noGenpart()) { // Generator level
230 #ifdef EDM_ML_DEBUG
231  edm::LogVerbatim("CherenkovAnalysis") << "XtalDedxAnalysis::Track found with ID " << simTrkItr->type();
232 #endif
233  mType_->Fill(simTrkItr->type());
234  }
235  break;
236  }
237  }
238  }
239  }
240  }
241 }
242 
243 // define this as a plug-in
XtalDedxAnalysis::caloHitSource_
edm::InputTag caloHitSource_
Definition: XtalDedxAnalysis.cc:52
XtalDedxAnalysis::tok_calo_
edm::EDGetTokenT< edm::PCaloHitContainer > tok_calo_
Definition: XtalDedxAnalysis.cc:55
mps_fire.i
i
Definition: mps_fire.py:355
MessageLogger.h
hfClusterShapes_cfi.hits
hits
Definition: hfClusterShapes_cfi.py:5
ESHandle.h
edm::EDGetTokenT< edm::PCaloHitContainer >
XtalDedxAnalysis::meE1T0_
TH1F * meE1T0_[4]
Definition: XtalDedxAnalysis.cc:59
edm::LogInfo
Definition: MessageLogger.h:254
XtalDedxAnalysis::simTkLabel_
std::string simTkLabel_
Definition: XtalDedxAnalysis.cc:53
XtalDedxAnalysis::~XtalDedxAnalysis
~XtalDedxAnalysis() override
Definition: XtalDedxAnalysis.cc:40
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
EDAnalyzer.h
XtalDedxAnalysis::meNHit_
TH1F * meNHit_[4]
Definition: XtalDedxAnalysis.cc:59
edm::Handle< edm::SimTrackContainer >
XtalDedxAnalysis::analyzeHits
void analyzeHits(std::vector< PCaloHit > &, edm::Handle< edm::SimTrackContainer > &, edm::Handle< edm::SimVertexContainer > &)
Definition: XtalDedxAnalysis.cc:144
edm::EDAnalyzer
Definition: EDAnalyzer.h:29
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Service.h
SimVertex.h
tfile
Definition: tfile.py:1
XtalDedxAnalysis::meE9T0_
TH1F * meE9T0_[4]
Definition: XtalDedxAnalysis.cc:59
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
XtalDedxAnalysis::endJob
void endJob() override
Definition: XtalDedxAnalysis.cc:45
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
TFileService.h
edm::ParameterSet
Definition: ParameterSet.h:36
Event.h
XtalDedxAnalysis::meE9T1_
TH1F * meE9T1_[4]
Definition: XtalDedxAnalysis.cc:59
PCaloHit.h
XtalDedxAnalysis::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: XtalDedxAnalysis.cc:119
edm::Service< TFileService >
edm::LogVerbatim
Definition: MessageLogger.h:297
edm::EventSetup
Definition: EventSetup.h:57
XtalDedxAnalysis::tok_vtx_
edm::EDGetTokenT< edm::SimVertexContainer > tok_vtx_
Definition: XtalDedxAnalysis.cc:57
XtalDedxAnalysis
Definition: XtalDedxAnalysis.cc:37
InputTag.h
compare.tfile
tfile
Definition: compare.py:325
overlapproblemtsosanalyzer_cfi.title
title
Definition: overlapproblemtsosanalyzer_cfi.py:7
types
types
Definition: AlignPCLThresholds_PayloadInspector.cc:30
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
XtalDedxAnalysis::beginJob
void beginJob() override
Definition: XtalDedxAnalysis.cc:43
type
type
Definition: HCALResponse.h:21
Frameworkfwd.h
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
EventSetup.h
Exception.h
PCaloHitContainer.h
cms::Exception
Definition: Exception.h:70
XtalDedxAnalysis::XtalDedxAnalysis
XtalDedxAnalysis(const edm::ParameterSet &)
Definition: XtalDedxAnalysis.cc:63
XtalDedxAnalysis::mType_
TH1I * mType_
Definition: XtalDedxAnalysis.cc:60
SimTrack.h
ParameterSet.h
XtalDedxAnalysis::meE1T1_
TH1F * meE1T1_[4]
Definition: XtalDedxAnalysis.cc:59
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
ntuplemaker.time
time
Definition: ntuplemaker.py:310
edm::Event
Definition: Event.h:73
XtalDedxAnalysis::tok_tk_
edm::EDGetTokenT< edm::SimTrackContainer > tok_tk_
Definition: XtalDedxAnalysis.cc:56
class-composition.parent
parent
Definition: class-composition.py:88
SimTrackContainer.h
edm::InputTag
Definition: InputTag.h:15
SimVertexContainer.h
hit
Definition: SiStripHitEffFromCalibTree.cc:88
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37