CMS 3D CMS Logo

EcalSimHitDump.cc
Go to the documentation of this file.
5 
15 
18 
19 #include <string>
20 #include <vector>
21 
23 public:
25  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
26 
27 protected:
28  void analyze(edm::Event const&, edm::EventSetup const&) override;
29 
30 private:
32  const std::vector<std::string> hitLab_;
33  const std::vector<edm::EDGetTokenT<edm::PCaloHitContainer>> toksCalo_;
34  const std::vector<int> types_;
35  const int maxEvent_;
36  int kount_;
37 };
38 
40  : g4Label_(ps.getParameter<std::string>("ModuleLabel")),
41  hitLab_(ps.getParameter<std::vector<std::string>>("HitCollections")),
42  toksCalo_{edm::vector_transform(
43  hitLab_,
44  [this](const std::string& name) { return consumes<edm::PCaloHitContainer>(edm::InputTag{g4Label_, name}); })},
45  types_(ps.getParameter<std::vector<int>>("CollectionTypes")),
46  maxEvent_(ps.getParameter<int>("MaxEvent")),
47  kount_(0) {
48  edm::LogVerbatim("HitStudy") << "Module Label: " << g4Label_ << " with " << hitLab_.size()
49  << " collections and maxEvent = " << maxEvent_;
50  for (unsigned int k = 0; k < hitLab_.size(); ++k)
51  edm::LogVerbatim("HitStudy") << "[" << k << "] Type " << types_[k] << " Label " << hitLab_[k];
52 }
53 
56  std::vector<std::string> coll = {"EcalHitsEB", "EcalHitsEE", "EcalHitsES"};
57  std::vector<int> type = {0, 1, 2};
58  desc.add<std::string>("ModuleLabel", "g4SimHits");
59  desc.add<std::vector<std::string>>("HitCollections", coll);
60  desc.add<std::vector<int>>("CollectionTypes", type);
61  desc.add<int>("MaxEvent", 10);
62  descriptions.add("ecalSimHitDump", desc);
63 }
64 
66  ++kount_;
67  edm::LogVerbatim("HitStudy") << "[" << kount_ << "] Run = " << e.id().run() << " Event = " << e.id().event();
68 
69  if ((kount_ <= maxEvent_) || (maxEvent_ <= 0)) {
70  for (unsigned int k = 0; k < toksCalo_.size(); ++k) {
71  const edm::Handle<edm::PCaloHitContainer>& hitsCalo = e.getHandle(toksCalo_[k]);
72  if (hitsCalo.isValid())
73  edm::LogVerbatim("HitStudy") << "EcalSimHitDump: Input " << hitsCalo->size() << " hits of type " << types_[k];
74  unsigned int i(0);
75  for (auto const& hit : *hitsCalo) {
76  double edep = hit.energy();
77  double time = hit.time();
78  unsigned int id = hit.id();
79  if (types_[k] == 0)
80  edm::LogVerbatim("HitStudy") << "[" << i << "] " << EBDetId(id) << " E" << edep << " T " << time;
81  else if (types_[k] == 1)
82  edm::LogVerbatim("HitStudy") << "[" << i << "] " << EEDetId(id) << " E" << edep << " T " << time;
83  else
84  edm::LogVerbatim("HitStudy") << "[" << i << "] " << ESDetId(id) << " E" << edep << " T " << time;
85  ++i;
86  }
87  }
88  }
89 }
90 
91 //define this as a plug-in
Log< level::Info, true > LogVerbatim
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const std::vector< int > types_
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
const int maxEvent_
EcalSimHitDump(const edm::ParameterSet &ps)
const std::vector< edm::EDGetTokenT< edm::PCaloHitContainer > > toksCalo_
const std::string g4Label_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
unsigned int id
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool isValid() const
Definition: HandleBase.h:70
void analyze(edm::Event const &, edm::EventSetup const &) override
const std::vector< std::string > hitLab_