CMS 3D CMS Logo

HGCalTimingAnalyzer.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 #include <iostream>
4 #include <fstream>
5 #include <vector>
6 #include <map>
7 #include <string>
8 
9 // user include files
16 
27 
40 
41 // Root objects
42 #include "TROOT.h"
43 #include "TSystem.h"
44 #include "TTree.h"
45 
46 //#define EDM_ML_DEBUG
47 
48 class HGCalTimingAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns,edm::one::SharedResources> {
49 
50 public:
51  explicit HGCalTimingAnalyzer(edm::ParameterSet const&);
53 
54  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
55 
56 private:
57  virtual void beginJob() override ;
58  virtual void beginRun(edm::Run const&, edm::EventSetup const&) override;
59  virtual void endRun(edm::Run const&, edm::EventSetup const&) override {}
60  virtual void analyze(edm::Event const&, edm::EventSetup const&) override;
61  void analyzeSimHits(int type, std::vector<PCaloHit> const& hits);
64 
69  double timeUnit_;
70  std::vector<int> idBeams_;
75  TTree *tree_;
76  std::vector<uint32_t> simHitCellIdEE_, simHitCellIdBeam_;
77  std::vector<float> simHitCellEnEE_, simHitCellEnBeam_;
78  std::vector<float> simHitCellTmEE_, simHitCellTmBeam_;
80 };
81 
83 
84  usesResource("TFileService");
85 
86  //now do whatever initialization is needed
87  detectorEE_ = iConfig.getParameter<std::string>("DetectorEE");
88  detectorBeam_= iConfig.getParameter<std::string>("DetectorBeam");
89  // Group hits (if groupHits_ = true) if hits come within timeUnit_
90  groupHits_ = iConfig.getParameter<bool>("GroupHits");
91  timeUnit_ = iConfig.getParameter<double>("TimeUnit");
92  // Only look into the beam counters with ID's as in idBeams_
93  idBeams_ = iConfig.getParameter<std::vector<int>>("IDBeams");
94  doTree_ = iConfig.getUntrackedParameter<bool>("DoTree",false);
95  if (!groupHits_) timeUnit_ = 0.000001;
96 #ifdef EDM_ML_DEBUG
97  std::cout << "HGCalTimingAnalyzer:: Group Hits " << groupHits_ << " in "
98  << timeUnit_ << " IdBeam " << idBeams_.size() << ":";
99  for (const auto& id : idBeams_) std::cout << " " << id;
100  std::cout << std::endl;
101 #endif
102  if (idBeams_.size() == 0) idBeams_.push_back(1001);
103 
104  edm::InputTag tmp0 = iConfig.getParameter<edm::InputTag>("GeneratorSrc");
105  tok_hepMC_ = consumes<edm::HepMCProduct>(tmp0);
106 #ifdef EDM_ML_DEBUG
107  std::cout << "HGCalTimingAnalyzer:: GeneratorSource = " << tmp0 << std::endl;
108 #endif
109  tok_simTk_ = consumes<edm::SimTrackContainer>(edm::InputTag("g4SimHits"));
110  tok_simVtx_ = consumes<edm::SimVertexContainer>(edm::InputTag("g4SimHits"));
111  std::string tmp1 = iConfig.getParameter<std::string>("CaloHitSrcEE");
112  tok_hitsEE_ = consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits",tmp1));
113 #ifdef EDM_ML_DEBUG
114  std::cout << "HGCalTimingAnalyzer:: Detector " << detectorEE_
115  << " with tags " << tmp1 << std::endl;
116 #endif
117  tmp1 = iConfig.getParameter<std::string>("CaloHitSrcBeam");
118  tok_hitsBeam_= consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits",tmp1));
119 #ifdef EDM_ML_DEBUG
120  std::cout << "HGCalTimingAnalyzer:: Detector " << detectorBeam_
121  << " with tags " << tmp1 << std::endl;
122 #endif
123 }
124 
126 
129  desc.add<std::string>("DetectorEE","HGCalEESensitive");
130  desc.add<std::string>("DetectorBeam","HcalTB06BeamDetector");
131  desc.add<bool>("GroupHits",false);
132  desc.add<double>("TimeUnit",0.001);
133  std::vector<int> ids = {1001,1002,1003,1004,1005};
134  desc.add<std::vector<int>>("IDBeams",ids);
135  desc.addUntracked<bool>("DoTree",true);
136  desc.add<edm::InputTag>("GeneratorSrc",edm::InputTag("generatorSmeared"));
137  desc.add<std::string>("CaloHitSrcEE","HGCHitsEE");
138  desc.add<std::string>("CaloHitSrcBeam","HcalTB06BeamHits");
139  descriptions.add("HGCalTimingAnalyzer",desc);
140 }
141 
143  std::string det(detectorEE_);
144  if (doTree_) {
145  tree_ = fs_->make<TTree>("HGCTB","SimHitEnergy");
146  tree_->Branch("xBeam", &xBeam_, "xBeam/D");
147  tree_->Branch("yBeam", &yBeam_, "yBeam/D");
148  tree_->Branch("zBeam", &zBeam_, "zBeam/D");
149  tree_->Branch("pBeam", &pBeam_, "pBeam/D");
150  tree_->Branch("simHitCellIdEE_", &simHitCellIdEE_);
151  tree_->Branch("simHitCellEnEE_", &simHitCellEnEE_);
152  tree_->Branch("simHitCellTmEE_", &simHitCellTmEE_);
153  tree_->Branch("simHitCellIdBeam_", &simHitCellIdBeam_);
154  tree_->Branch("simHitCellEnBeam_", &simHitCellEnBeam_);
155  tree_->Branch("simHitCellTmBeam_", &simHitCellTmBeam_);
156  }
157 }
158 
161  iSetup.get<IdealGeometryRecord>().get(detectorEE_, pHGDC);
162  hgcons_ = &(*pHGDC);
163 
164 #ifdef EDM_ML_DEBUG
165  std::cout << "HGCalTimingAnalyzer::" << detectorEE_ << " defined with "
166  << hgcons_->layers(false) << " layers" << std::endl;
167 #endif
168 
169 }
170 
172  const edm::EventSetup& iSetup) {
173 
174 #ifdef EDM_ML_DEBUG
175  //Generator input
177  iEvent.getByToken(tok_hepMC_,evtMC);
178  if (!evtMC.isValid()) {
179  edm::LogWarning("HGCal") << "no HepMCProduct found";
180  } else {
181  const HepMC::GenEvent * myGenEvent = evtMC->GetEvent();
182  unsigned int k(0);
183  for (HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
184  p != myGenEvent->particles_end(); ++p, ++k) {
185  std::cout << "Particle[" << k << "] with p " << (*p)->momentum().rho()
186  << " theta " << (*p)->momentum().theta() << " phi "
187  << (*p)->momentum().phi() << std::endl;
188  }
189  }
190 #endif
191 
192  //Now the Simhits
194  iEvent.getByToken(tok_simTk_, SimTk);
196  iEvent.getByToken(tok_simVtx_, SimVtx);
197  analyzeSimTracks(SimTk, SimVtx);
198 
199  simHitCellIdEE_.clear(); simHitCellIdBeam_.clear();
200  simHitCellEnEE_.clear(); simHitCellEnBeam_.clear();
201  simHitCellTmEE_.clear(); simHitCellTmBeam_.clear();
202 
203  edm::Handle<edm::PCaloHitContainer> theCaloHitContainers;
204  std::vector<PCaloHit> caloHits;
205  iEvent.getByToken(tok_hitsEE_, theCaloHitContainers);
206  if (theCaloHitContainers.isValid()) {
207 #ifdef EDM_ML_DEBUG
208  std::cout << "PcalohitContainer for " << detectorEE_ << " has "
209  << theCaloHitContainers->size() << " hits" << std::endl;
210 #endif
211  caloHits.clear();
212  caloHits.insert(caloHits.end(), theCaloHitContainers->begin(),
213  theCaloHitContainers->end());
214  analyzeSimHits(0, caloHits);
215  } else {
216 #ifdef EDM_ML_DEBUG
217  std::cout << "PCaloHitContainer does not exist for " << detectorEE_
218  << " !!!" << std::endl;
219 #endif
220  }
221 
222  iEvent.getByToken(tok_hitsBeam_, theCaloHitContainers);
223  if (theCaloHitContainers.isValid()) {
224 #ifdef EDM_ML_DEBUG
225  std::cout << "PcalohitContainer for " << detectorBeam_ << " has "
226  << theCaloHitContainers->size() << " hits" << std::endl;
227 #endif
228  caloHits.clear();
229  caloHits.insert(caloHits.end(), theCaloHitContainers->begin(),
230  theCaloHitContainers->end());
231  analyzeSimHits(1, caloHits);
232  } else {
233 #ifdef EDM_ML_DEBUG
234  std::cout << "PCaloHitContainer does not exist for " << detectorBeam_
235  << " !!!" << std::endl;
236 #endif
237  }
238  if (doTree_) tree_->Fill();
239 
240 }
241 
243  std::vector<PCaloHit> const& hits) {
244 
245 #ifdef EDM_ML_DEBUG
246  unsigned int i(0);
247 #endif
248  std::map<std::pair<uint32_t,uint64_t>,std::pair<double,double> > map_hits;
249  for (const auto& hit : hits) {
250  double energy = hit.energy();
251  double time = hit.time();
252  uint32_t id = hit.id();
253  if (type == 0) {
254  int subdet, zside, layer, sector, subsector, cell;
255  HGCalTestNumbering::unpackHexagonIndex(id, subdet, zside, layer, sector,
256  subsector, cell);
257  std::pair<int,int> recoLayerCell = hgcons_->simToReco(cell,layer,sector,true);
258  id = HGCalDetId((ForwardSubdetector)(subdet),zside,
259  recoLayerCell.second,subsector,sector,
260  recoLayerCell.first).rawId();
261 #ifdef EDM_ML_DEBUG
262  std::cout << "SimHit:Hit[" << i << "] Id " << subdet << ":" << zside
263  << ":" << layer << ":" << sector << ":" << subsector << ":"
264  << recoLayerCell.first << ":" << recoLayerCell.second
265  << " Energy " << energy << " Time " << time << std::endl;
266 #endif
267  } else {
268 #ifdef EDM_ML_DEBUG
269  int subdet, layer, x, y;
270  HcalTestBeamNumbering::unpackIndex(id, subdet, layer, x, y);
271  std::cout << "SimHit:Hit[" << i << "] Beam Subdet " << subdet
272  << " Layer " << layer << " x|y " << x << ":" << y
273  << " Energy " << energy << " Time " << time << std::endl;
274 #endif
275  }
276  uint64_t tid = (uint64_t)((time+50.0)/timeUnit_);
277  std::pair<uint32_t,uint64_t> key(id,tid);
278  auto itr = map_hits.find(key);
279  if (itr == map_hits.end()) {
280  map_hits[key] = std::pair<double,double>(time,0.0);
281  itr = map_hits.find(key);
282  }
283  energy += (itr->second).second;
284  map_hits[key] = std::pair<double,double>((itr->second).first,energy);
285 #ifdef EDM_ML_DEBUG
286  ++i;
287 #endif
288  }
289 
290 #ifdef EDM_ML_DEBUG
291  std::cout << "analyzeSimHits: Finds " << map_hits.size() << " hits "
292  << " from the Hit Vector of size " << hits.size() << " for type "
293  << type << std::endl;
294 #endif
295  for (const auto& itr: map_hits) {
296  uint32_t id = (itr.first).first;
297  double time = (itr.second).first;
298  double energy = (itr.second).second;
299  if (type == 0) {
300  simHitCellIdEE_.push_back(id);
301  simHitCellEnEE_.push_back(energy);
302  simHitCellTmEE_.push_back(time);
303  } else {
304  simHitCellIdBeam_.push_back(id);
305  simHitCellEnBeam_.push_back(energy);
306  simHitCellTmBeam_.push_back(time);
307  }
308 #ifdef EDM_ML_DEBUG
309  std::cout << "SimHit::ID: " << std::hex << id << std::dec << " T: " << time
310  << " E: " << energy << std::endl;
311 #endif
312  }
313 }
314 
315 
317  edm::Handle<edm::SimVertexContainer> const& SimVtx) {
318 
319  xBeam_ = yBeam_ = zBeam_ = pBeam_ = -1000000;
320  int vertIndex(-1);
321  for (edm::SimTrackContainer::const_iterator simTrkItr = SimTk->begin();
322  simTrkItr!= SimTk->end(); simTrkItr++) {
323 #ifdef EDM_ML_DEBUG
324  std::cout << "Track " << simTrkItr->trackId() << " Vertex "
325  << simTrkItr->vertIndex() << " Type " << simTrkItr->type()
326  << " Charge " << simTrkItr->charge() << " momentum "
327  << simTrkItr->momentum() << " " << simTrkItr->momentum().P()
328  << std::endl;
329 #endif
330  if (vertIndex == -1) {
331  vertIndex = simTrkItr->vertIndex();
332  pBeam_ = simTrkItr->momentum().P();
333  }
334  }
335  if (vertIndex != -1 && vertIndex < (int)SimVtx->size()) {
336  edm::SimVertexContainer::const_iterator simVtxItr= SimVtx->begin();
337  for (int iv=0; iv<vertIndex; iv++) simVtxItr++;
338 #ifdef EDM_ML_DEBUG
339  std::cout << "Vertex " << vertIndex << " position "
340  << simVtxItr->position() << std::endl;
341 #endif
342  xBeam_ = simVtxItr->position().X();
343  yBeam_ = simVtxItr->position().Y();
344  zBeam_ = simVtxItr->position().Z();
345  }
346 
347 }
348 
349 //define this as a plug-in
350 
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
const HGCalDDDConstants * hgcons_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
virtual void analyze(edm::Event const &, edm::EventSetup const &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< edm::SimTrackContainer > tok_simTk_
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hitsEE_
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
ForwardSubdetector
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
U second(std::pair< T, U > const &p)
HGCalTimingAnalyzer(edm::ParameterSet const &)
edm::EDGetTokenT< edm::HepMCProduct > tok_hepMC_
std::pair< int, int > simToReco(int cell, int layer, int mod, bool half) const
int iEvent
Definition: GenABIO.cc:230
std::vector< float > simHitCellEnEE_
unsigned int layers(bool reco) const
std::vector< float > simHitCellTmBeam_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
std::vector< int > idBeams_
std::vector< float > simHitCellEnBeam_
std::vector< uint32_t > simHitCellIdBeam_
std::vector< float > simHitCellTmEE_
std::vector< uint32_t > simHitCellIdEE_
virtual void beginJob() override
int k[5][pyjets_maxn]
unsigned int id
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
edm::Service< TFileService > fs_
unsigned long long uint64_t
Definition: Time.h:15
edm::EDGetTokenT< edm::SimVertexContainer > tok_simVtx_
const T & get() const
Definition: EventSetup.h:55
void add(std::string const &label, ParameterSetDescription const &psetDescription)
virtual void beginRun(edm::Run const &, edm::EventSetup const &) override
static void unpackIndex(const uint32_t &idx, int &det, int &lay, int &x, int &y)
void analyzeSimHits(int type, std::vector< PCaloHit > const &hits)
virtual void endRun(edm::Run const &, edm::EventSetup const &) override
static void unpackHexagonIndex(const uint32_t &idx, int &subdet, int &z, int &lay, int &wafer, int &celltyp, int &cell)
void analyzeSimTracks(edm::Handle< edm::SimTrackContainer > const &SimTk, edm::Handle< edm::SimVertexContainer > const &SimVtx)
Definition: Run.h:42
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hitsBeam_