1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
|
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "Geometry/CaloGeometry/interface/CaloGeometry.h"
#include "Geometry/Records/interface/CaloGeometryRecord.h"
#include "Geometry/HcalTowerAlgo/interface/HcalGeometry.h"
#include <iomanip>
#include <iostream>
class HcalCellParameterDump : public edm::one::EDAnalyzer<> {
public:
explicit HcalCellParameterDump(const edm::ParameterSet&);
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
void beginJob() override {}
void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
void endJob() override {}
private:
static constexpr int detMax_ = 4;
int subdet_;
edm::ESGetToken<CaloGeometry, CaloGeometryRecord> tok_geom_;
};
HcalCellParameterDump::HcalCellParameterDump(const edm::ParameterSet& iConfig) {
subdet_ = std::min(detMax_, std::max(iConfig.getParameter<int>("SubDetector"), 1));
subdet_ = std::min(detMax_, std::max(subdet_, 1));
tok_geom_ = esConsumes<CaloGeometry, CaloGeometryRecord>();
edm::LogVerbatim("HCalGeom") << "Dumps all cells for SubDetId: " << subdet_;
}
void HcalCellParameterDump::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<int>("SubDetector", 1);
descriptions.add("hcalCellParameterDump", desc);
}
void HcalCellParameterDump::analyze(const edm::Event& /*iEvent*/, const edm::EventSetup& iSetup) {
const CaloGeometry* geo = &iSetup.getData(tok_geom_);
const HcalGeometry* hcalGeom = static_cast<const HcalGeometry*>(geo->getSubdetectorGeometry(DetId::Hcal, HcalBarrel));
std::string subdets[detMax_] = {"HB", "HE", "HO", "HF"};
HcalSubdetector subdetd[detMax_] = {HcalBarrel, HcalEndcap, HcalOuter, HcalForward};
edm::LogVerbatim("HCalGeom") << "\n\nStudy Detector = Hcal SubDetector = " << subdets[subdet_ - 1]
<< "\n======================================\n";
const std::vector<DetId>& ids = hcalGeom->getValidDetIds(DetId::Hcal, subdetd[subdet_ - 1]);
int nall(0);
for (auto id : ids) {
++nall;
auto geom = hcalGeom->getGeometry(id);
edm::LogVerbatim("HCalGeom") << "[" << nall << "] " << HcalDetId(id) << " Reference " << std::setprecision(4)
<< geom->getPosition() << " Back " << geom->getBackPoint() << " [r,eta,phi] ("
<< geom->rhoPos() << ", " << geom->etaPos() << ":" << geom->etaSpan() << ", "
<< geom->phiPos() << ":" << geom->phiSpan() << ")";
}
edm::LogVerbatim("HCalGeom") << "\n\nDumps " << nall << " cells of the detector\n";
}
DEFINE_FWK_MODULE(HcalCellParameterDump);
|