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