CMS 3D CMS Logo

HcalGeomCheck.cc
Go to the documentation of this file.
1 // system include files
2 #include <cmath>
3 #include <iostream>
4 #include <fstream>
5 #include <vector>
6 #include <map>
7 #include <string>
8 
10 
12 
22 
27 
30 
31 #include "CLHEP/Geometry/Point3D.h"
32 #include "CLHEP/Geometry/Transform3D.h"
33 #include "CLHEP/Geometry/Vector3D.h"
34 #include "CLHEP/Units/GlobalSystemOfUnits.h"
35 #include "CLHEP/Units/GlobalPhysicalConstants.h"
36 
37 #include "TH1D.h"
38 #include "TH2D.h"
39 
40 class HcalGeomCheck : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
41 public:
42  struct hitsinfo {
43  hitsinfo() { phi = eta = energy = time = 0.0; }
44  double phi, eta, energy, time;
45  };
46 
47  explicit HcalGeomCheck(const edm::ParameterSet&);
48  ~HcalGeomCheck() override = default;
49  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
50 
51 protected:
52  void beginJob() override;
53  void beginRun(edm::Run const&, edm::EventSetup const&) override;
54  void analyze(edm::Event const&, edm::EventSetup const&) override;
55  void endRun(edm::Run const&, edm::EventSetup const&) override {}
56 
57 private:
58  void analyzeHits(int, const std::string&, const std::vector<PCaloHit>&);
59 
60  // ----------member data ---------------------------
63  const double rmin_, rmax_, zmin_, zmax_;
64  const int nbinR_, nbinZ_, verbosity_;
68 
69  //histogram related stuff
70  TH2D* h_RZ_;
71  std::map<std::pair<int, int>, TH1D*> h_phi_;
72  std::map<int, TH1D*> h_eta_;
73  std::vector<TH1D*> h_E_, h_T_;
74 };
75 
77  : caloHitSource_(iConfig.getParameter<std::string>("caloHitSource")),
78  ietaMin_(iConfig.getUntrackedParameter<int>("ietaMin", -41)),
79  ietaMax_(iConfig.getUntrackedParameter<int>("ietaMax", 41)),
80  depthMax_(iConfig.getUntrackedParameter<int>("depthMax", 7)),
81  rmin_(iConfig.getUntrackedParameter<double>("rMin", 0.0)),
82  rmax_(iConfig.getUntrackedParameter<double>("rMax", 5500.0)),
83  zmin_(iConfig.getUntrackedParameter<double>("zMin", -12500.0)),
84  zmax_(iConfig.getUntrackedParameter<double>("zMax", 12500.0)),
85  nbinR_(iConfig.getUntrackedParameter<int>("nBinR", 550)),
86  nbinZ_(iConfig.getUntrackedParameter<int>("nBinZ", 2500)),
87  verbosity_(iConfig.getUntrackedParameter<int>("verbosity", 0)),
88  tok_hits_(consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits", caloHitSource_))),
90  usesResource(TFileService::kSharedResource);
91 }
92 
95  desc.add<std::string>("caloHitSource", "HcalHits");
96  desc.addUntracked<int>("ietaMin", -41);
97  desc.addUntracked<int>("ietaMax", 41);
98  desc.addUntracked<int>("depthMax", 7);
99  desc.addUntracked<double>("rMin", 0.0);
100  desc.addUntracked<double>("rMax", 5500.0);
101  desc.addUntracked<double>("zMin", -12500.0);
102  desc.addUntracked<double>("zMax", 12500.0);
103  desc.addUntracked<int>("nBinR", 550);
104  desc.addUntracked<int>("nBinZ", 250);
105  desc.addUntracked<int>("verbosity", 0);
106  descriptions.add("hcalGeomCheck", desc);
107 }
108 
110  //Now the hits
111  const edm::Handle<edm::PCaloHitContainer>& theCaloHitContainer = iEvent.getHandle(tok_hits_);
112  if (theCaloHitContainer.isValid()) {
113  if (verbosity_ > 0)
114  edm::LogVerbatim("HcalValidation") << " PcalohitItr = " << theCaloHitContainer->size();
115 
116  //Merge hits for the same DetID
117  std::map<HcalDetId, hitsinfo> map_hits;
118  unsigned int nused(0);
119  for (auto const& hit : *(theCaloHitContainer.product())) {
120  unsigned int id = hit.id();
122  double energy = hit.energy();
123  double time = hit.time();
124  int subdet = detId.subdet();
125  int ieta = detId.ieta();
126  int iphi = detId.iphi();
127  int depth = detId.depth();
128  std::pair<double, double> etaphi = hcons_->getEtaPhi(subdet, ieta, iphi);
129  double rz = hcons_->getRZ(subdet, ieta, depth);
130  if (verbosity_ > 2)
131  edm::LogVerbatim("HcalValidation") << "i/p " << subdet << ":" << ieta << ":" << iphi << ":" << depth << " o/p "
132  << etaphi.first << ":" << etaphi.second << ":" << rz;
133  HepGeom::Point3D<float> gcoord = HepGeom::Point3D<float>(rz * cos(etaphi.second) / cosh(etaphi.first),
134  rz * sin(etaphi.second) / cosh(etaphi.first),
135  rz * tanh(etaphi.first));
136  nused++;
137  double tof = (gcoord.mag() * CLHEP::mm) / CLHEP::c_light;
138  if (verbosity_ > 1)
139  edm::LogVerbatim("HcalValidation")
140  << "Detector " << subdet << " ieta = " << ieta << " iphi = " << iphi << " depth = " << depth
141  << " positon = " << gcoord << " energy = " << energy << " time = " << time << ":" << tof;
142  time -= tof;
143  if (time < 0)
144  time = 0;
145  hitsinfo hinfo;
146  if (map_hits.count(detId) != 0) {
147  hinfo = map_hits[detId];
148  } else {
149  hinfo.phi = gcoord.getPhi();
150  hinfo.eta = gcoord.getEta();
151  hinfo.time = time;
152  }
153  hinfo.energy += energy;
154  map_hits[detId] = hinfo;
155 
156  h_RZ_->Fill(gcoord.z(), gcoord.rho());
157  }
158 
159  //Fill in histograms
160  for (auto const& hit : map_hits) {
161  hitsinfo hinfo = hit.second;
162  int subdet = hit.first.subdet();
163  if (subdet > 0 && subdet < static_cast<int>(h_E_.size())) {
164  int depth = hit.first.depth();
165  int ieta = hit.first.ieta();
166  int iphi = hit.first.iphi();
167  if (verbosity_ > 1)
168  edm::LogVerbatim("HGCalValidation")
169  << " ---------------------- eta = " << ieta << ":" << hinfo.eta << " phi = " << iphi << ":" << hinfo.phi
170  << " depth = " << depth << " E = " << hinfo.energy << " T = " << hinfo.time;
171  h_E_[subdet]->Fill(hinfo.energy);
172  h_T_[subdet]->Fill(hinfo.time);
173  auto itr1 = h_phi_.find(std::pair<int, int>(ieta, depth));
174  if (itr1 != h_phi_.end())
175  itr1->second->Fill(iphi);
176  auto itr2 = h_eta_.find(depth);
177  if (itr2 != h_eta_.end())
178  itr2->second->Fill(ieta);
179  }
180  }
181  } else if (verbosity_ > 0) {
182  edm::LogVerbatim("HcalValidation") << "PCaloHitContainer does not "
183  << "exist for HCAL";
184  }
185 }
186 
187 // ------------ method called when starting to processes a run ------------
188 void HcalGeomCheck::beginRun(const edm::Run&, const edm::EventSetup& iSetup) {
189  const auto& pHRNDC = iSetup.getData(tok_HRNDC_);
190  hcons_ = &pHRNDC;
191  if (verbosity_ > 0)
192  edm::LogVerbatim("HcalValidation") << "Obtain HcalDDDRecConstants from Event Setup";
193 }
194 
196  if (verbosity_ > 2)
197  edm::LogVerbatim("HcalValidation") << "HcalGeomCheck:: Enter beginJob";
199  h_RZ_ = fs->make<TH2D>("RZ", "R vs Z", nbinZ_, zmin_, zmax_, nbinR_, rmin_, rmax_);
200  if (verbosity_ > 2)
201  edm::LogVerbatim("HcalValidation") << "HcalGeomCheck: booked scatterplot RZ";
202  char name[20], title[100];
203  for (int depth = 1; depth < depthMax_; ++depth) {
204  for (int ieta = ietaMin_; ieta <= ietaMax_; ++ieta) {
205  sprintf(name, "phi%d%d", ieta, depth);
206  sprintf(title, "i#phi (i#eta = %d, depth = %d)", ieta, depth);
207  h_phi_[std::pair<int, int>(ieta, depth)] = fs->make<TH1D>(name, title, 400, -20, 380);
208  if (verbosity_ > 2)
209  edm::LogVerbatim("HcalValidation") << "HcalGeomCheck: books " << title;
210  }
211  sprintf(name, "eta%d", depth);
212  sprintf(title, "i#eta (depth = %d)", depth);
213  h_eta_[depth] = fs->make<TH1D>(name, title, 100, -50, 50);
214  if (verbosity_ > 2)
215  edm::LogVerbatim("HcalValidation") << "HcalGeomCheck: books " << title;
216  }
217 
218  std::vector<std::string> dets = {"HCAL", "HB", "HE", "HF", "HO"};
219  for (unsigned int ih = 0; ih < dets.size(); ++ih) {
220  sprintf(name, "E_%s", dets[ih].c_str());
221  sprintf(title, "Energy deposit in %s (MeV)", dets[ih].c_str());
222  h_E_.emplace_back(fs->make<TH1D>(name, title, 1000, 0.0, 1.0));
223  if (verbosity_ > 2)
224  edm::LogVerbatim("HcalValidation") << "HcalGeomCheck: books " << title;
225  sprintf(name, "T_%s", dets[ih].c_str());
226  sprintf(title, "Time of hit in %s (ns)", dets[ih].c_str());
227  h_T_.emplace_back(fs->make<TH1D>(name, title, 1000, 0.0, 200.0));
228  if (verbosity_ > 2)
229  edm::LogVerbatim("HcalValidation") << "HcalGeomCheck: books " << title;
230  }
231 }
232 
234 //define this as a plug-in
std::pair< T, T > etaphi(T x, T y, T z)
Definition: FastMath.h:162
static const std::string kSharedResource
Definition: TFileService.h:76
Log< level::Info, true > LogVerbatim
const double rmax_
double getRZ(const int &subdet, const int &ieta, const int &depth) const
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
std::vector< PCaloHit > PCaloHitContainer
std::vector< TH1D * > h_T_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void tanh(data_T data[CONFIG_T::n_in], res_T res[CONFIG_T::n_in])
void beginJob() override
const int ietaMin_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T const * product() const
Definition: Handle.h:70
void beginRun(edm::Run const &, edm::EventSetup const &) override
~HcalGeomCheck() override=default
const HcalDDDRecConstants * hcons_
const int depthMax_
const int nbinR_
std::pair< double, double > getEtaPhi(const int &subdet, const int &ieta, const int &iphi) const
void analyze(edm::Event const &, edm::EventSetup const &) override
constexpr HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:138
int iEvent
Definition: GenABIO.cc:224
const int ietaMax_
std::map< std::pair< int, int >, TH1D * > h_phi_
constexpr int ieta() const
get the cell ieta
Definition: HcalDetId.h:155
void analyzeHits(int, const std::string &, const std::vector< PCaloHit > &)
std::map< int, TH1D * > h_eta_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
const edm::ESGetToken< HcalDDDRecConstants, HcalRecNumberingRecord > tok_HRNDC_
Transition
Definition: Transition.h:12
const double zmin_
std::vector< TH1D * > h_E_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool getData(T &iHolder) const
Definition: EventSetup.h:122
void endRun(edm::Run const &, edm::EventSetup const &) override
DetId relabel(const uint32_t testId) const
unsigned int id
const int nbinZ_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< edm::PCaloHitContainer > tok_hits_
const std::string caloHitSource_
bool isValid() const
Definition: HandleBase.h:70
HLT enums.
const double rmin_
const double zmax_
HcalGeomCheck(const edm::ParameterSet &)
constexpr int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
Definition: Run.h:45
const int verbosity_
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164