CMS 3D CMS Logo

L1CaloInputScaleTester.cc
Go to the documentation of this file.
2 
3 // system include files
4 #include <memory>
5 #include <iostream>
6 using std::cout;
7 using std::endl;
8 #include <iomanip>
9 using std::setprecision;
10 
11 // user include files
13 
17 
19 
28 
29 //
30 // constants, enums and typedefs
31 //
32 
33 //
34 // static data member definitions
35 //
36 
37 //
38 // constructors and destructor
39 //
41  : hcalScaleToken_(esConsumes<L1CaloHcalScale, L1CaloHcalScaleRcd>()),
42  ecalScaleToken_(esConsumes<L1CaloEcalScale, L1CaloEcalScaleRcd>()),
43  transcoderToken_(esConsumes<CaloTPGTranscoder, CaloTPGRecord>()),
44  tokens_(consumesCollector()) {
45  //now do what ever initialization is needed
46 }
47 
49  // do anything here that needs to be done at destruction time
50  // (e.g. close files, deallocate resources etc.)
51 }
52 
53 //
54 // member functions
55 //
56 
57 // ------------ method called to for each event ------------
59  using namespace edm;
60 
61 #ifdef THIS_IS_AN_EVENT_EXAMPLE
63  iEvent.getByLabel("example", pIn);
64 #endif
65 
66 #ifdef THIS_IS_AN_EVENTSETUP_EXAMPLE
67  ESHandle<SetupData> pSetup;
68  iSetup.get<SetupRecord>().get(pSetup);
69 #endif
70 
71  ESHandle<L1CaloEcalScale> caloEcalScale = iSetup.getHandle(ecalScaleToken_);
72  ESHandle<L1CaloHcalScale> caloHcalScale = iSetup.getHandle(hcalScaleToken_);
73  ESHandle<CaloTPGTranscoder> caloTPGTranscoder = iSetup.getHandle(transcoderToken_);
74 
75  EcalTPGScale ecalTPGScale(tokens_, iSetup);
76 
77  bool ecalIsConsistent = true;
78  bool hcalIsConsistent = true;
79 
80  double ecal1;
81  double ecal2;
82  double hcal1;
83  double hcal2;
84  double hcal3;
85  double hcal4;
86 
87  // compare the ecal scales
88 
89  // 8 bits of input energy
90  for (unsigned short input = 0; input <= 0xFF; input++) {
91  // loop over ietas, barrel first
92  for (unsigned short absIeta = 1; absIeta <= 17; absIeta++) {
93  // positive eta
94  ecal1 = ecalTPGScale.getTPGInGeV((unsigned int)input, EcalTrigTowerDetId(1, EcalBarrel, absIeta, 1));
95  ecal2 = caloEcalScale->et(input, absIeta, 1);
96 
97  if (!(ecal1 == ecal2)) {
98  ecalIsConsistent = false;
99  /*LogWarning("InconsistentData")
100  << "ECAL scales not consistent! (pos eta barrel)"
101  << "old ECAL is " << setprecision (8) << ecal1
102  << " new ECAL is " << setprecision (8) << ecal2
103  << " difference is " << (ecal1 - ecal2) ; */
104  }
105 
106  // negative eta
107  ecal1 = ecalTPGScale.getTPGInGeV((unsigned int)input, EcalTrigTowerDetId(-1, EcalBarrel, absIeta, 2));
108  ecal2 = caloEcalScale->et(input, absIeta, -1);
109 
110  if (!(ecal1 == ecal2)) {
111  ecalIsConsistent = false;
112  /*LogWarning("InconsistentData")
113  << "ECAL scales not consistent! (neg eta barrel)"
114  << "old ECAL is " << setprecision (8) << ecal1
115  << " new ECAL is " << setprecision (8) << ecal2
116  << " difference is " << (ecal1 - ecal2) ; */
117  }
118  }
119  // now loop over endcap ietas
120  for (unsigned short absIeta = 18; absIeta <= 28; absIeta++) {
121  // positive eta
122  ecal1 = ecalTPGScale.getTPGInGeV((unsigned int)input, EcalTrigTowerDetId(1, EcalEndcap, absIeta, 1));
123  ecal2 = caloEcalScale->et(input, absIeta, 1);
124 
125  if (!(ecal1 == ecal2)) {
126  ecalIsConsistent = false;
127  /*LogWarning("InconsistentData")
128  << "ECAL scales not consistent! (pos eta endcap)"
129  << "old ECAL is " << setprecision (8) << ecal1
130  << " new ECAL is " << setprecision (8) << ecal2
131  << " difference is " << (ecal1 - ecal2) ; */
132  }
133 
134  // negative eta
135  ecal1 = ecalTPGScale.getTPGInGeV((unsigned int)input, EcalTrigTowerDetId(-1, EcalEndcap, absIeta, 2));
136  ecal2 = caloEcalScale->et(input, absIeta, -1);
137 
138  if (!(ecal1 == ecal2)) {
139  ecalIsConsistent = false;
140  /*LogWarning("InconsistentData")
141  << "ECAL scales not consistent! (neg eta endcap)"
142  << "old ECAL is " << setprecision (8) << ecal1
143  << " new ECAL is " << setprecision (8) << ecal2
144  << " difference is " << (ecal1 - ecal2) ; */
145  }
146  }
147  }
148 
149  if (!ecalIsConsistent) {
150  // do something
151  //cout << "WARNING: ECAL scales not consistent!" << endl;
152  LogWarning("InconsistentData") << "ECAL scales not consistent!";
153  } else {
154  // do something else
155  //cout << "ECAL scales okay" << endl;
156  }
157 
158  // compare the hcal scales
159 
160  for (unsigned short input = 0; input <= 0xFF; input++) {
161  // loop over ietas
162  for (unsigned short absIeta = 1; absIeta <= 32; absIeta++) {
163  hcal1 = caloTPGTranscoder->hcaletValue(absIeta, input); // no eta-
164  hcal2 = caloHcalScale->et(input, absIeta, 1); // sign in transcoder
165  hcal3 = caloTPGTranscoder->hcaletValue(-absIeta, input); // no eta-
166  hcal4 = caloHcalScale->et(input, absIeta, -1); // sign in transcoder
167 
168  if ((!(hcal1 == hcal2)) || (!(hcal3 == hcal4))) {
169  hcalIsConsistent = false;
170  /*LogWarning("InconsistentData")
171  << "HCAL scales not consistent!"
172  << "old HCAL is " << hcal1
173  << " new HCAL is " << hcal2 ; */
174  }
175  }
176  }
177  if (!hcalIsConsistent) {
178  // do something
179  //cout << "WARNING: HCAL scales not consistent!" << endl;
180  LogWarning("InconsistentData") << "HCAL scales not consistent!";
181  } else {
182  // do something else
183  //cout << "HCAL scales okay" << endl;
184  }
185 }
186 
187 // ------------ method called once each job just before starting event loop ------------
189 
190 // ------------ method called once each job just after ending the event loop ------------
L1CaloInputScaleTester(const edm::ParameterSet &)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
double et(unsigned short rank, unsigned short eta, short etaSign) const
convert from rank to physically meaningful quantity
static std::string const input
Definition: EdmProvDump.cc:50
EcalTPGScale::Tokens tokens_
int iEvent
Definition: GenABIO.cc:224
T get() const
Definition: EventSetup.h:79
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
std::vector< edm::EDGetTokenT< int > > tokens_
double getTPGInGeV(const EcalTriggerPrimitiveDigi &tpDigi) const
Definition: EcalTPGScale.cc:17
HLT enums.
edm::ESGetToken< CaloTPGTranscoder, CaloTPGRecord > transcoderToken_
double et(unsigned short rank, unsigned short eta, short etaSign) const
convert from rank to physically meaningful quantity
Log< level::Warning, false > LogWarning
edm::ESGetToken< L1CaloEcalScale, L1CaloEcalScaleRcd > ecalScaleToken_
void analyze(const edm::Event &, const edm::EventSetup &) override
edm::ESGetToken< L1CaloHcalScale, L1CaloHcalScaleRcd > hcalScaleToken_