CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HGCalTimingAnalyzer.cc
Go to the documentation of this file.
1 // system include files
2 #include <fstream>
3 #include <sstream>
4 #include <map>
5 #include <string>
6 #include <vector>
7 
8 // user include files
13 
24 
36 
37 // Root objects
38 #include "TROOT.h"
39 #include "TSystem.h"
40 #include "TTree.h"
41 
42 //#define EDM_ML_DEBUG
43 
44 class HGCalTimingAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
45 public:
46  explicit HGCalTimingAnalyzer(edm::ParameterSet const&);
47  ~HGCalTimingAnalyzer() override;
48 
49  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
50 
51 private:
52  void beginJob() override;
53  void beginRun(edm::Run const&, edm::EventSetup const&) override;
54  void endRun(edm::Run const&, edm::EventSetup const&) override {}
55  void analyze(edm::Event const&, edm::EventSetup const&) override;
56  void analyzeSimHits(int type, std::vector<PCaloHit> const& hits);
59 
65  double timeUnit_;
66  std::vector<int> idBeams_;
71  TTree* tree_;
72  std::vector<uint32_t> simHitCellIdEE_, simHitCellIdBeam_;
73  std::vector<float> simHitCellEnEE_, simHitCellEnBeam_;
74  std::vector<float> simHitCellTmEE_, simHitCellTmBeam_;
76 };
77 
79  : detectorEE_(iConfig.getParameter<std::string>("DetectorEE")),
80  detectorBeam_(iConfig.getParameter<std::string>("DetectorBeam")),
82  edm::ESInputTag("", detectorEE_))) {
83  usesResource("TFileService");
84 
85  // now do whatever initialization is needed
86  // Group hits (if groupHits_ = true) if hits come within timeUnit_
87  groupHits_ = iConfig.getParameter<bool>("GroupHits");
88  timeUnit_ = iConfig.getParameter<double>("TimeUnit");
89  // Only look into the beam counters with ID's as in idBeams_
90  idBeams_ = iConfig.getParameter<std::vector<int>>("IDBeams");
91  doTree_ = iConfig.getUntrackedParameter<bool>("DoTree", false);
92  if (!groupHits_)
93  timeUnit_ = 0.000001;
94 #ifdef EDM_ML_DEBUG
95  std::ostringstream st1;
96  st1 << "HGCalTimingAnalyzer:: Group Hits " << groupHits_ << " in " << timeUnit_ << " IdBeam " << idBeams_.size()
97  << ":";
98  for (const auto& id : idBeams_)
99  st1 << " " << id;
100  edm::LogVerbatim("HGCSim") << st1.str();
101 #endif
102  if (idBeams_.empty())
103  idBeams_.push_back(1001);
104 
105  edm::InputTag tmp0 = iConfig.getParameter<edm::InputTag>("GeneratorSrc");
106  tok_hepMC_ = consumes<edm::HepMCProduct>(tmp0);
107 #ifdef EDM_ML_DEBUG
108  edm::LogVerbatim("HGCSim") << "HGCalTimingAnalyzer:: GeneratorSource = " << tmp0;
109 #endif
110  tok_simTk_ = consumes<edm::SimTrackContainer>(edm::InputTag("g4SimHits"));
111  tok_simVtx_ = consumes<edm::SimVertexContainer>(edm::InputTag("g4SimHits"));
112  std::string tmp1 = iConfig.getParameter<std::string>("CaloHitSrcEE");
113  tok_hitsEE_ = consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits", tmp1));
114 #ifdef EDM_ML_DEBUG
115  edm::LogVerbatim("HGCSim") << "HGCalTimingAnalyzer:: Detector " << detectorEE_ << " with tags " << tmp1;
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  edm::LogVerbatim("HGCSim") << "HGCalTimingAnalyzer:: Detector " << detectorBeam_ << " with tags " << tmp1;
121 #endif
122 }
123 
125 
128  desc.add<std::string>("DetectorEE", "HGCalEESensitive");
129  desc.add<std::string>("DetectorBeam", "HcalTB06BeamDetector");
130  desc.add<bool>("GroupHits", false);
131  desc.add<double>("TimeUnit", 0.001);
132  std::vector<int> ids = {1001, 1002, 1003, 1004, 1005};
133  desc.add<std::vector<int>>("IDBeams", ids);
134  desc.addUntracked<bool>("DoTree", true);
135  desc.add<edm::InputTag>("GeneratorSrc", edm::InputTag("generatorSmeared"));
136  desc.add<std::string>("CaloHitSrcEE", "HGCHitsEE");
137  desc.add<std::string>("CaloHitSrcBeam", "HcalTB06BeamHits");
138  descriptions.add("HGCalTimingAnalyzer", desc);
139 }
140 
143  if (doTree_) {
144  tree_ = fs_->make<TTree>("HGCTB", "SimHitEnergy");
145  tree_->Branch("xBeam", &xBeam_, "xBeam/D");
146  tree_->Branch("yBeam", &yBeam_, "yBeam/D");
147  tree_->Branch("zBeam", &zBeam_, "zBeam/D");
148  tree_->Branch("pBeam", &pBeam_, "pBeam/D");
149  tree_->Branch("simHitCellIdEE_", &simHitCellIdEE_);
150  tree_->Branch("simHitCellEnEE_", &simHitCellEnEE_);
151  tree_->Branch("simHitCellTmEE_", &simHitCellTmEE_);
152  tree_->Branch("simHitCellIdBeam_", &simHitCellIdBeam_);
153  tree_->Branch("simHitCellEnBeam_", &simHitCellEnBeam_);
154  tree_->Branch("simHitCellTmBeam_", &simHitCellTmBeam_);
155  }
156 }
157 
159  hgcons_ = &iSetup.getData(tokDDD_);
160 
161 #ifdef EDM_ML_DEBUG
162  edm::LogVerbatim("HGCSim") << "HGCalTimingAnalyzer::" << detectorEE_ << " defined with " << hgcons_->layers(false)
163  << " layers";
164 #endif
165 }
166 
168 #ifdef EDM_ML_DEBUG
169  // Generator input
171  iEvent.getByToken(tok_hepMC_, evtMC);
172  if (!evtMC.isValid()) {
173  edm::LogWarning("HGCal") << "no HepMCProduct found";
174  } else {
175  const HepMC::GenEvent* myGenEvent = evtMC->GetEvent();
176  unsigned int k(0);
177  for (HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end();
178  ++p, ++k) {
179  edm::LogVerbatim("HGCSim") << "Particle[" << k << "] with p " << (*p)->momentum().rho() << " theta "
180  << (*p)->momentum().theta() << " phi " << (*p)->momentum().phi();
181  }
182  }
183 #endif
184 
185  // Now the Simhits
187  iEvent.getByToken(tok_simTk_, SimTk);
189  iEvent.getByToken(tok_simVtx_, SimVtx);
190  analyzeSimTracks(SimTk, SimVtx);
191 
192  simHitCellIdEE_.clear();
193  simHitCellIdBeam_.clear();
194  simHitCellEnEE_.clear();
195  simHitCellEnBeam_.clear();
196  simHitCellTmEE_.clear();
197  simHitCellTmBeam_.clear();
198 
199  edm::Handle<edm::PCaloHitContainer> theCaloHitContainers;
200  std::vector<PCaloHit> caloHits;
201  iEvent.getByToken(tok_hitsEE_, theCaloHitContainers);
202  if (theCaloHitContainers.isValid()) {
203 #ifdef EDM_ML_DEBUG
204  edm::LogVerbatim("HGCSim") << "PcalohitContainer for " << detectorEE_ << " has " << theCaloHitContainers->size()
205  << " hits";
206 #endif
207  caloHits.clear();
208  caloHits.insert(caloHits.end(), theCaloHitContainers->begin(), theCaloHitContainers->end());
209  analyzeSimHits(0, caloHits);
210  } else {
211 #ifdef EDM_ML_DEBUG
212  edm::LogVerbatim("HGCSim") << "PCaloHitContainer does not exist for " << detectorEE_ << " !!!";
213 #endif
214  }
215 
216  iEvent.getByToken(tok_hitsBeam_, theCaloHitContainers);
217  if (theCaloHitContainers.isValid()) {
218 #ifdef EDM_ML_DEBUG
219  edm::LogVerbatim("HGCSim") << "PcalohitContainer for " << detectorBeam_ << " has " << theCaloHitContainers->size()
220  << " hits";
221 #endif
222  caloHits.clear();
223  caloHits.insert(caloHits.end(), theCaloHitContainers->begin(), theCaloHitContainers->end());
224  analyzeSimHits(1, caloHits);
225  } else {
226 #ifdef EDM_ML_DEBUG
227  edm::LogVerbatim("HGCSim") << "PCaloHitContainer does not exist for " << detectorBeam_ << " !!!";
228 #endif
229  }
230  if (doTree_)
231  tree_->Fill();
232 }
233 
234 void HGCalTimingAnalyzer::analyzeSimHits(int type, std::vector<PCaloHit> const& hits) {
235 #ifdef EDM_ML_DEBUG
236  unsigned int i(0);
237 #endif
238  std::map<std::pair<uint32_t, uint64_t>, std::pair<double, double>> map_hits;
239  for (const auto& hit : hits) {
240  double energy = hit.energy();
241  double time = hit.time();
242  uint32_t id = hit.id();
243  if (type == 0) {
244  int subdet, zside, layer, sector, subsector, cell;
245  HGCalTestNumbering::unpackHexagonIndex(id, subdet, zside, layer, sector, subsector, cell);
246  std::pair<int, int> recoLayerCell = hgcons_->simToReco(cell, layer, sector, true);
247  id = HGCalDetId((ForwardSubdetector)(subdet), zside, recoLayerCell.second, subsector, sector, recoLayerCell.first)
248  .rawId();
249 #ifdef EDM_ML_DEBUG
250  edm::LogVerbatim("HGCSim") << "SimHit:Hit[" << i << "] Id " << subdet << ":" << zside << ":" << layer << ":"
251  << sector << ":" << subsector << ":" << recoLayerCell.first << ":"
252  << recoLayerCell.second << " Energy " << energy << " Time " << time;
253 #endif
254  } else {
255 #ifdef EDM_ML_DEBUG
256  int subdet, layer, x, y;
257  HcalTestBeamNumbering::unpackIndex(id, subdet, layer, x, y);
258  edm::LogVerbatim("HGCSim") << "SimHit:Hit[" << i << "] Beam Subdet " << subdet << " Layer " << layer << " x|y "
259  << x << ":" << y << " Energy " << energy << " Time " << time;
260 #endif
261  }
262  uint64_t tid = (uint64_t)((time + 50.0) / timeUnit_);
263  std::pair<uint32_t, uint64_t> key(id, tid);
264  auto itr = map_hits.find(key);
265  if (itr == map_hits.end()) {
266  map_hits[key] = std::pair<double, double>(time, 0.0);
267  itr = map_hits.find(key);
268  }
269  energy += (itr->second).second;
270  map_hits[key] = std::pair<double, double>((itr->second).first, energy);
271 #ifdef EDM_ML_DEBUG
272  ++i;
273 #endif
274  }
275 
276 #ifdef EDM_ML_DEBUG
277  edm::LogVerbatim("HGCSim") << "analyzeSimHits: Finds " << map_hits.size() << " hits "
278  << " from the Hit Vector of size " << hits.size() << " for type " << type;
279 #endif
280  for (const auto& itr : map_hits) {
281  uint32_t id = (itr.first).first;
282  double time = (itr.second).first;
283  double energy = (itr.second).second;
284  if (type == 0) {
285  simHitCellIdEE_.push_back(id);
286  simHitCellEnEE_.push_back(energy);
287  simHitCellTmEE_.push_back(time);
288  } else {
289  simHitCellIdBeam_.push_back(id);
290  simHitCellEnBeam_.push_back(energy);
291  simHitCellTmBeam_.push_back(time);
292  }
293 #ifdef EDM_ML_DEBUG
294  edm::LogVerbatim("HGCSim") << "SimHit::ID: " << std::hex << id << std::dec << " T: " << time << " E: " << energy;
295 #endif
296  }
297 }
298 
300  edm::Handle<edm::SimVertexContainer> const& SimVtx) {
301  xBeam_ = yBeam_ = zBeam_ = pBeam_ = -1000000;
302  int vertIndex(-1);
303  for (edm::SimTrackContainer::const_iterator simTrkItr = SimTk->begin(); simTrkItr != SimTk->end(); simTrkItr++) {
304 #ifdef EDM_ML_DEBUG
305  edm::LogVerbatim("HGCSim") << "Track " << simTrkItr->trackId() << " Vertex " << simTrkItr->vertIndex() << " Type "
306  << simTrkItr->type() << " Charge " << simTrkItr->charge() << " momentum "
307  << simTrkItr->momentum() << " " << simTrkItr->momentum().P();
308 #endif
309  if (vertIndex == -1) {
310  vertIndex = simTrkItr->vertIndex();
311  pBeam_ = simTrkItr->momentum().P();
312  }
313  }
314  if (vertIndex != -1 && vertIndex < (int)SimVtx->size()) {
315  edm::SimVertexContainer::const_iterator simVtxItr = SimVtx->begin();
316  for (int iv = 0; iv < vertIndex; iv++)
317  simVtxItr++;
318 #ifdef EDM_ML_DEBUG
319  edm::LogVerbatim("HGCSim") << "Vertex " << vertIndex << " position " << simVtxItr->position();
320 #endif
321  xBeam_ = simVtxItr->position().X();
322  yBeam_ = simVtxItr->position().Y();
323  zBeam_ = simVtxItr->position().Z();
324  }
325 }
326 
327 // define this as a plug-in
328 
Log< level::Info, true > LogVerbatim
T getUntrackedParameter(std::string const &, T const &) const
int32_t *__restrict__ iv
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
const HGCalDDDConstants * hgcons_
uint16_t *__restrict__ id
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
void analyze(edm::Event const &, edm::EventSetup const &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::ESGetToken< HGCalDDDConstants, IdealGeometryRecord > tokDDD_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
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)
constexpr std::array< uint8_t, layerIndexSize > layer
bool getData(T &iHolder) const
Definition: EventSetup.h:128
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:224
std::vector< float > simHitCellEnEE_
unsigned int layers(bool reco) const
tuple key
prepare the HTCondor submission files and eventually submit them
Transition
Definition: Transition.h:12
std::vector< float > simHitCellTmBeam_
const std::string detectorBeam_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
std::vector< int > idBeams_
std::vector< float > simHitCellEnBeam_
std::vector< uint32_t > simHitCellIdBeam_
std::vector< float > simHitCellTmEE_
std::vector< uint32_t > simHitCellIdEE_
unsigned int id
void endRun(edm::Run const &, edm::EventSetup const &) override
edm::Service< TFileService > fs_
unsigned long long uint64_t
Definition: Time.h:13
edm::EDGetTokenT< edm::SimVertexContainer > tok_simVtx_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const std::string detectorEE_
void beginRun(edm::Run const &, edm::EventSetup const &) override
static void unpackIndex(const uint32_t &idx, int &det, int &lay, int &x, int &y)
Log< level::Warning, false > LogWarning
void analyzeSimHits(int type, std::vector< PCaloHit > const &hits)
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
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:45
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hitsBeam_