CMS 3D CMS Logo

HcalTestNumberingTest.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HcalTestNumberingTester
4 // Class: HcalTestNumberingTester
5 //
13 //
14 // Original Author: Sunanda Banerjee
15 // Created: Mon 2013/12/26
16 //
17 
18 // system include files
19 #include <memory>
20 #include <iostream>
21 #include <fstream>
22 
23 // user include files
28 
33 
44 
46 public:
48  ~HcalTestNumberingTester() override;
49 
50  void beginJob() override {}
51  void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
52  void endJob() override {}
53 };
54 
56 
58 
59 // ------------ method called to produce the data ------------
62  iSetup.get<HcalSimNumberingRecord>().get(pHSNDC);
64  iSetup.get<HcalRecNumberingRecord>().get(pHRNDC);
65 
66  if (pHSNDC.isValid() && pHRNDC.isValid()) {
67  std::cout << "about to de-reference the edm's" << std::endl;
68  HcalDDDSimConstants* hcs = (HcalDDDSimConstants*)(&(*pHSNDC));
69  HcalDDDRecConstants* hcr = (HcalDDDRecConstants*)(&(*pHRNDC));
71  HcalNumberingScheme* schme2 = dynamic_cast<HcalNumberingScheme*>(new HcalTestNumberingScheme(false));
72 
73  for (int type = 0; type < 2; ++type) {
74  HcalSubdetector sub = (type == 0) ? HcalBarrel : HcalEndcap;
75  for (int zs = 0; zs < 2; ++zs) {
76  int zside = 2 * zs - 1;
77  std::pair<int, int> etas = hcr->getEtaRange(type);
78  for (int eta = etas.first; eta <= etas.second; ++eta) {
79  std::vector<std::pair<int, double> > phis = hcr->getPhis(sub, eta);
80  for (unsigned int k = 0; k < phis.size(); ++k) {
81  int phi = phis[k].first;
82  int lmin = (type == 1 && eta == 16) ? 8 : 1;
83  int lmax = (type == 1) ? 19 : ((eta == 16) ? 7 : 17);
84  for (int lay = lmin; lay <= lmax; ++lay) {
85  std::pair<int, int> etd = hcs->getEtaDepth(sub, eta, phi, zside, 0, lay);
86  HcalNumberingFromDDD::HcalID tmp(sub, zs, etd.second, etd.first, phi, phi, lay);
87  uint32_t id1 = schme1->getUnitID(tmp);
88  uint32_t id2 = schme2->getUnitID(tmp);
90  std::cout << "I/P " << sub << ":" << zside * eta << ":" << phi << ":" << lay << " Normal " << std::hex
91  << id1 << std::dec << " " << HcalDetId(id1) << " Test " << std::hex << id2 << std::dec << " "
92  << HcalDetId(id0);
93  if (id1 != id0.rawId())
94  std::cout << " *** ERROR ***";
95  std::cout << std::endl;
96  }
97  }
98  }
99  }
100  }
101  }
102 }
103 
104 //define this as a plug-in
EDAnalyzer.h
HcalNumberingScheme::getUnitID
virtual uint32_t getUnitID(const HcalNumberingFromDDD::HcalID &id)
Definition: HcalNumberingScheme.cc:19
HcalDDDRecConstants::getPhis
std::vector< std::pair< int, double > > getPhis(const int &subdet, const int &ieta) const
Definition: HcalDDDRecConstants.cc:373
HcalNumberingFromDDD::HcalID
Definition: HcalNumberingFromDDD.h:21
ESHandle.h
HcalSimNumberingRecord.h
ecaldqm::zside
int zside(DetId const &)
Definition: EcalDQMCommonUtils.cc:189
globals_cff.id1
id1
Definition: globals_cff.py:33
gather_cfg.cout
cout
Definition: gather_cfg.py:144
HcalRecNumberingRecord.h
edm::one::EDAnalyzer
Definition: EDAnalyzer.h:30
HcalTestNumberingTester::~HcalTestNumberingTester
~HcalTestNumberingTester() override
Definition: HcalTestNumberingTest.cc:57
HcalBarrel
Definition: HcalAssistant.h:33
createJobs.tmp
tmp
align.sh
Definition: createJobs.py:716
HcalNumberingFromDDD.h
HcalDDDSimConstants
Definition: HcalDDDSimConstants.h:24
DetId
Definition: DetId.h:17
HcalTestNumberingTester::analyze
void analyze(edm::Event const &iEvent, edm::EventSetup const &) override
Definition: HcalTestNumberingTest.cc:60
MakerMacros.h
edm::EventSetup::get
T get() const
Definition: EventSetup.h:80
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
EnergyCorrector.etas
etas
Definition: EnergyCorrector.py:45
PVValHelper::eta
Definition: PVValidationHelpers.h:69
edm::ESHandle
Definition: DTSurvey.h:22
HcalNumberingScheme.h
dqmdumpme.k
k
Definition: dqmdumpme.py:60
HcalRecNumberingRecord
Definition: HcalRecNumberingRecord.h:23
HcalDDDSimConstants::getEtaDepth
std::pair< int, int > getEtaDepth(const int &det, int etaR, const int &phi, const int &zside, int depth, const int &lay) const
Definition: HcalDDDSimConstants.cc:252
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
HcalTestNumberingScheme
Definition: HcalTestNumberingScheme.h:11
type
type
Definition: SiPixelVCal_PayloadInspector.cc:37
HcalDetId.h
gainCalibHelper::gainCalibPI::type
type
Definition: SiPixelGainCalibHelper.h:39
HcalHitRelabeller.h
HcalDetId
Definition: HcalDetId.h:12
iEvent
int iEvent
Definition: GenABIO.cc:224
HcalTestNumberingTester
Definition: HcalTestNumberingTest.cc:45
IdealGeometryRecord.h
HcalTestNumberingTester::endJob
void endJob() override
Definition: HcalTestNumberingTest.cc:52
edm::EventSetup
Definition: EventSetup.h:57
HcalSubdetector.h
get
#define get
edm::ESHandleBase::isValid
bool isValid() const
Definition: ESHandle.h:44
HcalNumberingScheme
Definition: HcalNumberingScheme.h:13
HcalSubdetector
HcalSubdetector
Definition: HcalAssistant.h:31
DDAxes::phi
HcalHitRelabeller::relabel
DetId relabel(const uint32_t testId) const
Definition: HcalHitRelabeller.cc:49
DetId::rawId
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
HcalDDDRecConstants::getEtaRange
std::pair< int, int > getEtaRange(const int &i) const
Definition: HcalDDDRecConstants.h:74
HcalTestNumberingTester::beginJob
void beginJob() override
Definition: HcalTestNumberingTest.cc:50
HcalEndcap
Definition: HcalAssistant.h:34
Frameworkfwd.h
HcalDDDRecConstants.h
EventSetup.h
HcalDDDRecConstants
Definition: HcalDDDRecConstants.h:23
ParameterSet.h
HcalTestNumberingTester::HcalTestNumberingTester
HcalTestNumberingTester(const edm::ParameterSet &)
Definition: HcalTestNumberingTest.cc:55
globals_cff.id2
id2
Definition: globals_cff.py:34
edm::Event
Definition: Event.h:73
HcalTestNumberingScheme.h
HcalSimNumberingRecord
Definition: HcalSimNumberingRecord.h:25
TauDecayModes.dec
dec
Definition: TauDecayModes.py:143
HcalDDDSimConstants.h