CMS 3D CMS Logo

HcalSimHitDump.cc
Go to the documentation of this file.
3 
7 
13 
21 
23 
24 #include <memory>
25 #include <iostream>
26 #include <fstream>
27 #include <vector>
28 #include <map>
29 #include <string>
30 
32 public:
34  ~HcalSimHitDump() override = default;
35  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
36 
37 protected:
38  void beginJob() override {}
39  void endJob() override {}
40  void analyze(const edm::Event& e, const edm::EventSetup& c) override;
41 
42  void analyzeHits(std::vector<PCaloHit>&);
43 
44 private:
46  const int maxEvent_;
47  const bool testNumber_;
49  int nevt_;
50 };
51 
53  : g4Label_(ps.getParameter<std::string>("ModuleLabel")),
54  hitLab_(ps.getParameter<std::string>("HCCollection")),
55  maxEvent_(ps.getParameter<int>("MaxEvent")),
56  testNumber_(ps.getParameter<bool>("TestNumber")),
57  toks_calo_(consumes<edm::PCaloHitContainer>(edm::InputTag(g4Label_, hitLab_))),
58  nevt_(0) {
59  edm::LogVerbatim("HitStudy") << "HcalSimHitDump::Module Label: " << g4Label_ << " Hits: " << hitLab_ << " MaxEvent "
60  << maxEvent_ << " TestNumbering " << testNumber_;
61 }
62 
65  desc.add<std::string>("ModuleLabel", "g4SimHits");
66  desc.add<std::string>("HCCollection", "HcalHits");
67  desc.add<int>("MaxEvent", 10);
68  desc.add<bool>("TestNumber", true);
69  descriptions.add("hcalSimHitDump", desc);
70 }
71 
73  ++nevt_;
74  edm::LogVerbatim("HitStudy") << "HcalSimHitDump::Serial # " << nevt_ << " Run # " << e.id().run() << " Event # "
75  << e.id().event();
76 
77  if (nevt_ <= maxEvent_) {
78  std::vector<PCaloHit> hcHits;
79  const edm::Handle<edm::PCaloHitContainer>& hitsCalo = e.getHandle(toks_calo_);
80  if (hitsCalo.isValid()) {
81  edm::LogVerbatim("HitStudy") << "HcalValidation: get valid hist for Hcal";
82  std::vector<PCaloHit> caloHits;
83  caloHits.insert(caloHits.end(), hitsCalo->begin(), hitsCalo->end());
84  edm::LogVerbatim("HitStudy") << "HcalValidation: Hit buffer " << caloHits.size();
85  analyzeHits(caloHits);
86  }
87  }
88 }
89 
90 void HcalSimHitDump::analyzeHits(std::vector<PCaloHit>& hits) {
91  //Now the dump
92  for (unsigned int i = 0; i < hits.size(); i++) {
93  double edep = hits[i].energy();
94  double time = hits[i].time();
95  unsigned int id_ = hits[i].id();
96  if (testNumber_) {
97  int det, z, depth, eta, phi, lay;
99  std::string sub("HX");
100  if (det == 1)
101  sub = "HB";
102  else if (det == 2)
103  sub = "HE";
104  else if (det == 3)
105  sub = "HO";
106  else if (det == 4)
107  sub = "HF";
108  else if (det == 5)
109  sub = "HT";
110  int side = (z == 0) ? (-1) : (1);
111  edm::LogVerbatim("HitStudy") << "[" << i << "] (" << sub << " " << side * eta << "," << phi << "," << depth << ","
112  << lay << ") E " << edep << " T " << time;
113  } else {
114  edm::LogVerbatim("HitStudy") << "[" << i << "] " << HcalDetId(id_) << " E " << edep << " T " << time;
115  }
116  }
117 }
118 
119 //define this as a plug-in
Log< level::Info, true > LogVerbatim
std::vector< PCaloHit > PCaloHitContainer
const std::string g4Label_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void endJob() override
const bool testNumber_
const std::string hitLab_
void beginJob() override
static void unpackHcalIndex(const uint32_t &idx, int &det, int &z, int &depth, int &eta, int &phi, int &lay)
void analyzeHits(std::vector< PCaloHit > &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
~HcalSimHitDump() override=default
const int maxEvent_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool isValid() const
Definition: HandleBase.h:70
HLT enums.
const edm::EDGetTokenT< edm::PCaloHitContainer > toks_calo_
void analyze(const edm::Event &e, const edm::EventSetup &c) override
HcalSimHitDump(const edm::ParameterSet &ps)