CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 {}
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  :
78 
79  caloHitSource_(iConfig.getParameter<std::string>("caloHitSource")),
80  ietaMin_(iConfig.getUntrackedParameter<int>("ietaMin", -41)),
81  ietaMax_(iConfig.getUntrackedParameter<int>("ietaMax", 41)),
82  depthMax_(iConfig.getUntrackedParameter<int>("depthMax", 7)),
83  rmin_(iConfig.getUntrackedParameter<double>("rMin", 0.0)),
84  rmax_(iConfig.getUntrackedParameter<double>("rMax", 5500.0)),
85  zmin_(iConfig.getUntrackedParameter<double>("zMin", -12500.0)),
86  zmax_(iConfig.getUntrackedParameter<double>("zMax", 12500.0)),
87  nbinR_(iConfig.getUntrackedParameter<int>("nBinR", 550)),
88  nbinZ_(iConfig.getUntrackedParameter<int>("nBinZ", 2500)),
89  verbosity_(iConfig.getUntrackedParameter<int>("verbosity", 0)) {
90  usesResource(TFileService::kSharedResource);
91  tok_hits_ = consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits", caloHitSource_));
92  tok_HRNDC_ = esConsumes<HcalDDDRecConstants, HcalRecNumberingRecord, edm::Transition::BeginRun>();
93 }
94 
97  desc.add<std::string>("caloHitSource", "HcalHits");
98  desc.addUntracked<int>("ietaMin", -41);
99  desc.addUntracked<int>("ietaMax", 41);
100  desc.addUntracked<int>("depthMax", 7);
101  desc.addUntracked<double>("rMin", 0.0);
102  desc.addUntracked<double>("rMax", 5500.0);
103  desc.addUntracked<double>("zMin", -12500.0);
104  desc.addUntracked<double>("zMax", 12500.0);
105  desc.addUntracked<int>("nBinR", 550);
106  desc.addUntracked<int>("nBinZ", 250);
107  desc.addUntracked<int>("verbosity", 0);
108  descriptions.add("hcalGeomCheck", desc);
109 }
110 
112  //Now the hits
113  edm::Handle<edm::PCaloHitContainer> theCaloHitContainer;
114  iEvent.getByToken(tok_hits_, theCaloHitContainer);
115  if (theCaloHitContainer.isValid()) {
116  if (verbosity_ > 0)
117  edm::LogVerbatim("HcalValidation") << " PcalohitItr = " << theCaloHitContainer->size();
118 
119  //Merge hits for the same DetID
120  std::map<HcalDetId, hitsinfo> map_hits;
121  unsigned int nused(0);
122  for (auto const& hit : *(theCaloHitContainer.product())) {
123  unsigned int id = hit.id();
125  double energy = hit.energy();
126  double time = hit.time();
127  int subdet = detId.subdet();
128  int ieta = detId.ieta();
129  int iphi = detId.iphi();
130  int depth = detId.depth();
131  std::pair<double, double> etaphi = hcons_->getEtaPhi(subdet, ieta, iphi);
132  double rz = hcons_->getRZ(subdet, ieta, depth);
133  if (verbosity_ > 2)
134  edm::LogVerbatim("HcalValidation") << "i/p " << subdet << ":" << ieta << ":" << iphi << ":" << depth << " o/p "
135  << etaphi.first << ":" << etaphi.second << ":" << rz;
136  HepGeom::Point3D<float> gcoord = HepGeom::Point3D<float>(rz * cos(etaphi.second) / cosh(etaphi.first),
137  rz * sin(etaphi.second) / cosh(etaphi.first),
138  rz * tanh(etaphi.first));
139  nused++;
140  double tof = (gcoord.mag() * CLHEP::mm) / CLHEP::c_light;
141  if (verbosity_ > 1)
142  edm::LogVerbatim("HcalValidation")
143  << "Detector " << subdet << " ieta = " << ieta << " iphi = " << iphi << " depth = " << depth
144  << " positon = " << gcoord << " energy = " << energy << " time = " << time << ":" << tof;
145  time -= tof;
146  if (time < 0)
147  time = 0;
148  hitsinfo hinfo;
149  if (map_hits.count(detId) != 0) {
150  hinfo = map_hits[detId];
151  } else {
152  hinfo.phi = gcoord.getPhi();
153  hinfo.eta = gcoord.getEta();
154  hinfo.time = time;
155  }
156  hinfo.energy += energy;
157  map_hits[detId] = hinfo;
158 
159  h_RZ_->Fill(gcoord.z(), gcoord.rho());
160  }
161 
162  //Fill in histograms
163  for (auto const& hit : map_hits) {
164  hitsinfo hinfo = hit.second;
165  int subdet = hit.first.subdet();
166  if (subdet > 0 && subdet < static_cast<int>(h_E_.size())) {
167  int depth = hit.first.depth();
168  int ieta = hit.first.ieta();
169  int iphi = hit.first.iphi();
170  if (verbosity_ > 1)
171  edm::LogVerbatim("HGCalValidation")
172  << " ---------------------- eta = " << ieta << ":" << hinfo.eta << " phi = " << iphi << ":" << hinfo.phi
173  << " depth = " << depth << " E = " << hinfo.energy << " T = " << hinfo.time;
174  h_E_[subdet]->Fill(hinfo.energy);
175  h_T_[subdet]->Fill(hinfo.time);
176  auto itr1 = h_phi_.find(std::pair<int, int>(ieta, depth));
177  if (itr1 != h_phi_.end())
178  itr1->second->Fill(iphi);
179  auto itr2 = h_eta_.find(depth);
180  if (itr2 != h_eta_.end())
181  itr2->second->Fill(ieta);
182  }
183  }
184  } else if (verbosity_ > 0) {
185  edm::LogVerbatim("HcalValidation") << "PCaloHitContainer does not "
186  << "exist for HCAL";
187  }
188 }
189 
190 // ------------ method called when starting to processes a run ------------
191 void HcalGeomCheck::beginRun(const edm::Run&, const edm::EventSetup& iSetup) {
192  const auto& pHRNDC = iSetup.getData(tok_HRNDC_);
193  hcons_ = &pHRNDC;
194  if (verbosity_ > 0)
195  edm::LogVerbatim("HcalValidation") << "Obtain HcalDDDRecConstants from Event Setup";
196 }
197 
199  if (verbosity_ > 2)
200  edm::LogVerbatim("HcalValidation") << "HcalGeomCheck:: Enter beginJob";
202  h_RZ_ = fs->make<TH2D>("RZ", "R vs Z", nbinZ_, zmin_, zmax_, nbinR_, rmin_, rmax_);
203  if (verbosity_ > 2)
204  edm::LogVerbatim("HcalValidation") << "HcalGeomCheck: booked scatterplot RZ";
205  char name[20], title[100];
206  for (int depth = 1; depth < depthMax_; ++depth) {
207  for (int ieta = ietaMin_; ieta <= ietaMax_; ++ieta) {
208  sprintf(name, "phi%d%d", ieta, depth);
209  sprintf(title, "i#phi (i#eta = %d, depth = %d)", ieta, depth);
210  h_phi_[std::pair<int, int>(ieta, depth)] = fs->make<TH1D>(name, title, 400, -20, 380);
211  if (verbosity_ > 2)
212  edm::LogVerbatim("HcalValidation") << "HcalGeomCheck: books " << title;
213  }
214  sprintf(name, "eta%d", depth);
215  sprintf(title, "i#eta (depth = %d)", depth);
216  h_eta_[depth] = fs->make<TH1D>(name, title, 100, -50, 50);
217  if (verbosity_ > 2)
218  edm::LogVerbatim("HcalValidation") << "HcalGeomCheck: books " << title;
219  }
220 
221  std::vector<std::string> dets = {"HCAL", "HB", "HE", "HF", "HO"};
222  for (unsigned int ih = 0; ih < dets.size(); ++ih) {
223  sprintf(name, "E_%s", dets[ih].c_str());
224  sprintf(title, "Energy deposit in %s (MeV)", dets[ih].c_str());
225  h_E_.emplace_back(fs->make<TH1D>(name, title, 1000, 0.0, 1.0));
226  if (verbosity_ > 2)
227  edm::LogVerbatim("HcalValidation") << "HcalGeomCheck: books " << title;
228  sprintf(name, "T_%s", dets[ih].c_str());
229  sprintf(title, "Time of hit in %s (ns)", dets[ih].c_str());
230  h_T_.emplace_back(fs->make<TH1D>(name, title, 1000, 0.0, 200.0));
231  if (verbosity_ > 2)
232  edm::LogVerbatim("HcalValidation") << "HcalGeomCheck: books " << title;
233  }
234 }
235 
237 //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_
std::vector< TH1D * > h_T_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::pair< double, double > getEtaPhi(const int &subdet, const int &ieta, const int &iphi) const
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
void beginJob() override
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hits_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const int ietaMin_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
void beginRun(edm::Run const &, edm::EventSetup const &) override
const HcalDDDRecConstants * hcons_
const int depthMax_
const int nbinR_
bool getData(T &iHolder) const
Definition: EventSetup.h:128
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_
constexpr int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
std::map< std::pair< int, int >, TH1D * > h_phi_
void analyzeHits(int, const std::string &, const std::vector< PCaloHit > &)
double getRZ(const int &subdet, const int &ieta, const int &depth) const
std::map< int, TH1D * > h_eta_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
constexpr int ieta() const
get the cell ieta
Definition: HcalDetId.h:155
const double zmin_
std::vector< TH1D * > h_E_
void endRun(edm::Run const &, edm::EventSetup const &) override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
unsigned int id
~HcalGeomCheck() override
T const * product() const
Definition: Handle.h:70
const int nbinZ_
edm::ESGetToken< HcalDDDRecConstants, HcalRecNumberingRecord > tok_HRNDC_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const std::string caloHitSource_
const double rmin_
DetId relabel(const uint32_t testId) const
constexpr int depth() const
get the tower depth
Definition: HcalDetId.h:164
const double zmax_
HcalGeomCheck(const edm::ParameterSet &)
Definition: Run.h:45
const int verbosity_