CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HcalSimHitAnalysis.cc
Go to the documentation of this file.
3 
7 
13 
16 
21 
25 
32 
33 #include <TH2F.h>
34 #include <cmath>
35 #include <map>
36 #include <memory>
37 #include <string>
38 #include <vector>
39 
40 class HcalSimHitAnalysis : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
41 public:
43  ~HcalSimHitAnalysis() override {}
44  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
45 
46 protected:
47  void beginJob() override;
48  void analyze(edm::Event const&, edm::EventSetup const&) override;
49  void beginRun(edm::Run const&, edm::EventSetup const&) override {}
50  void endRun(edm::Run const&, edm::EventSetup const&) override {}
51 
52 private:
53  static const int ndets_ = 4;
60 };
61 
63  usesResource(TFileService::kSharedResource);
64 
65  g4Label_ = ps.getUntrackedParameter<std::string>("moduleLabel", "g4SimHits");
66  hcalHits_ = ps.getUntrackedParameter<std::string>("HitCollection", "HcalHits");
67  verbose_ = ps.getUntrackedParameter<bool>("Verbose", false);
68  testNumber_ = ps.getUntrackedParameter<bool>("TestNumber", true);
69 
70  tok_calo_ = consumes<edm::PCaloHitContainer>(edm::InputTag(g4Label_, hcalHits_));
71  tok_geom_ = esConsumes<CaloGeometry, CaloGeometryRecord>();
72  tok_HRNDC_ = esConsumes<HcalDDDRecConstants, HcalRecNumberingRecord>();
73 
74  edm::LogVerbatim("HitStudy") << "Module Label: " << g4Label_ << " Hits: " << hcalHits_ << " testNumber "
75  << testNumber_;
76 }
77 
80  desc.addUntracked<std::string>("ModuleLabel", "g4SimHits");
81  desc.addUntracked<std::string>("HitCollection", "HcalHits");
82  desc.addUntracked<bool>("Verbose", false);
83  desc.addUntracked<bool>("TestNumber", true);
84  descriptions.add("hcalSimHitAnalysis", desc);
85 }
86 
89  if (!tfile.isAvailable())
90  throw cms::Exception("BadConfig") << "TFileService unavailable: "
91  << "please add it to config file";
92  char name[20], title[120];
93  std::string dets[ndets_] = {"HB", "HE", "HO", "HF"};
94  int nx[ndets_] = {100, 100, 350, 160};
95  double xlo[ndets_] = {0, -300, 0, -160};
96  double xhi[ndets_] = {500, 300, 3500, 160};
97  int ny[ndets_] = {100, 100, 50, 160};
98  double ylo[ndets_] = {170, -300, 375, -160};
99  double yhi[ndets_] = {370, 300, 425, 160};
100  std::string xttl[ndets_] = {"|z| (cm)", "x (cm)", "|z| (cm)", "x (cm)"};
101  std::string yttl[ndets_] = {"#rho (cm)", "y (cm)", "#rho (cm)", "y (cm)"};
102  for (int i = 0; i < ndets_; i++) {
103  sprintf(name, "poszp%d", i);
104  sprintf(title, "%s+", dets[i].c_str());
105  poszp_[i] = tfile->make<TH2F>(name, title, nx[i], xlo[i], xhi[i], ny[i], ylo[i], yhi[i]);
106  poszp_[i]->GetXaxis()->SetTitle(xttl[i].c_str());
107  poszp_[i]->GetYaxis()->SetTitle(yttl[i].c_str());
108  sprintf(title, "%s-", dets[i].c_str());
109  poszp_[i]->GetYaxis()->SetTitleOffset(1.2);
110  poszp_[i]->Sumw2();
111  sprintf(name, "poszn%d", i);
112  poszn_[i] = tfile->make<TH2F>(name, title, nx[i], xlo[i], xhi[i], ny[i], ylo[i], yhi[i]);
113  poszn_[i]->GetXaxis()->SetTitle(xttl[i].c_str());
114  poszn_[i]->GetYaxis()->SetTitle(yttl[i].c_str());
115  poszn_[i]->GetYaxis()->SetTitleOffset(1.2);
116  poszn_[i]->Sumw2();
117  }
118 }
119 
121  if (verbose_)
122  edm::LogVerbatim("HitStudy") << "Run = " << e.id().run() << " Event = " << e.id().event();
123 
124  // get hcalGeometry
125  const CaloGeometry* geo = &iS.getData(tok_geom_);
126  const HcalGeometry* hgeom = static_cast<const HcalGeometry*>(geo->getSubdetectorGeometry(DetId::Hcal, HcalBarrel));
127  const auto& pHRNDC = iS.getData(tok_HRNDC_);
128  const HcalDDDRecConstants* hcons = &pHRNDC;
129 
131  e.getByToken(tok_calo_, hitsCalo);
132  bool getHits = (hitsCalo.isValid());
133  uint32_t nhits = (getHits) ? hitsCalo->size() : 0;
134  if (verbose_)
135  edm::LogVerbatim("HitStudy") << "HcalSimHitAnalysis: Input flags Hits " << getHits << " with " << nhits << " hits";
136  if (getHits) {
137  std::vector<PCaloHit> hits;
138  hits.insert(hits.end(), hitsCalo->begin(), hitsCalo->end());
139  if (!hits.empty()) {
140  std::map<HcalDetId, double> hitMap;
141  for (auto hit : hits) {
142  double edep = hit.energy();
143  uint32_t id = hit.id();
144  HcalDetId hid = (testNumber_) ? HcalHitRelabeller::relabel(id, hcons) : HcalDetId(id);
145  auto it = hitMap.find(hid);
146  if (it == hitMap.end()) {
147  hitMap[hid] = edep;
148  } else {
149  (it->second) += edep;
150  }
151  }
152 
153  for (auto it : hitMap) {
154  HcalDetId id(it.first);
155  GlobalPoint gpos = hgeom->getPosition(id);
156  HcalSubdetector subdet = (id).subdet();
157  int indx =
158  ((subdet == HcalBarrel)
159  ? 0
160  : ((subdet == HcalEndcap) ? 1 : ((subdet == HcalOuter) ? 2 : ((subdet == HcalForward) ? 3 : -1))));
161  if (verbose_)
162  edm::LogVerbatim("HitStudy") << "HcalSimHitAnalysis: " << id << ":" << it.second << " at " << gpos
163  << " subdet " << subdet << ":" << indx;
164  if (indx >= 0) {
165  double x = ((indx == 0) || (indx == 2)) ? std::abs(gpos.z()) : gpos.x();
166  double y = ((indx == 0) || (indx == 2)) ? (gpos.perp()) : gpos.y();
167  if (id.zside() >= 0)
168  poszp_[indx]->Fill(x, y);
169  else
170  poszn_[indx]->Fill(x, y);
171  }
172  }
173  }
174  }
175 }
176 
177 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:38
static const std::string kSharedResource
Definition: TFileService.h:76
Log< level::Info, true > LogVerbatim
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
EventNumber_t event() const
Definition: EventID.h:40
T getUntrackedParameter(std::string const &, T const &) const
HcalSimHitAnalysis(const edm::ParameterSet &ps)
void endRun(edm::Run const &, edm::EventSetup const &) override
T perp() const
Definition: PV3DBase.h:69
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
uint16_t *__restrict__ id
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T y() const
Definition: PV3DBase.h:60
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
int zside(DetId const &)
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > tok_geom_
bool getData(T &iHolder) const
Definition: EventSetup.h:128
void beginRun(edm::Run const &, edm::EventSetup const &) override
edm::EDGetTokenT< edm::PCaloHitContainer > tok_calo_
edm::ESGetToken< HcalDDDRecConstants, HcalRecNumberingRecord > tok_HRNDC_
T z() const
Definition: PV3DBase.h:61
bool isAvailable() const
Definition: Service.h:40
HcalSubdetector
Definition: HcalAssistant.h:31
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:70
static const int ndets_
GlobalPoint getPosition(const DetId &id) const
unsigned int id
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void analyze(edm::Event const &, edm::EventSetup const &) override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EventID id() const
Definition: EventBase.h:59
void beginJob() override
tuple tfile
Definition: compare.py:324
DetId relabel(const uint32_t testId) const
T x() const
Definition: PV3DBase.h:59
Definition: Run.h:45
~HcalSimHitAnalysis() override