CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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")),
43  [this](const std::string& name) {
44  return consumes<edm::PCaloHitContainer>(edm::InputTag{g4Label_, name});
45  })},
46  types_(ps.getParameter<std::vector<int>>("CollectionTypes")),
47  maxEvent_(ps.getParameter<int>("MaxEvent")),
48  kount_(0) {
49  edm::LogVerbatim("HitStudy") << "Module Label: " << g4Label_ << " with " << hitLab_.size()
50  << " collections and maxEvent = " << maxEvent_;
51  for (unsigned int k = 0; k < hitLab_.size(); ++k)
52  edm::LogVerbatim("HitStudy") << "[" << k << "] Type " << types_[k] << " Label " << hitLab_[k];
53 }
54 
57  std::vector<std::string> coll = {"EcalHitsEB", "EcalHitsEE", "EcalHitsES"};
58  std::vector<int> type = {0, 1, 2};
59  desc.add<std::string>("ModuleLabel", "g4SimHits");
60  desc.add<std::vector<std::string>>("HitCollections", coll);
61  desc.add<std::vector<int>>("CollectionTypes", type);
62  desc.add<int>("MaxEvent", 10);
63  descriptions.add("ecalSimHitDump", desc);
64 }
65 
67  ++kount_;
68  edm::LogVerbatim("HitStudy") << "[" << kount_ << "] Run = " << e.id().run() << " Event = " << e.id().event();
69 
70  if ((kount_ <= maxEvent_) || (maxEvent_ <= 0)) {
71  for (unsigned int k = 0; k < toksCalo_.size(); ++k) {
73  if (hitsCalo.isValid())
74  edm::LogVerbatim("HitStudy") << "EcalSimHitDump: Input " << hitsCalo->size() << " hits of type " << types_[k];
75  unsigned int i(0);
76  for (auto const& hit : *hitsCalo) {
77  double edep = hit.energy();
78  double time = hit.time();
79  unsigned int id = hit.id();
80  if (types_[k] == 0)
81  edm::LogVerbatim("HitStudy") << "[" << i << "] " << EBDetId(id) << " E" << edep << " T " << time;
82  else if (types_[k] == 1)
83  edm::LogVerbatim("HitStudy") << "[" << i << "] " << EEDetId(id) << " E" << edep << " T " << time;
84  else
85  edm::LogVerbatim("HitStudy") << "[" << i << "] " << ESDetId(id) << " E" << edep << " T " << time;
86  ++i;
87  }
88  }
89  }
90 }
91 
92 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:38
Log< level::Info, true > LogVerbatim
EventNumber_t event() const
Definition: EventID.h:40
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:563
const int maxEvent_
EcalSimHitDump(const edm::ParameterSet &ps)
const std::vector< edm::EDGetTokenT< edm::PCaloHitContainer > > toksCalo_
const std::string g4Label_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
unsigned int id
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EventID id() const
Definition: EventBase.h:59
void analyze(edm::Event const &, edm::EventSetup const &) override
const std::vector< std::string > hitLab_