CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 //
31 // constants, enums and typedefs
32 //
33 
34 //
35 // static data member definitions
36 //
37 
38 //
39 // constructors and destructor
40 //
42 
43 {
44  //now do what ever initialization is needed
45 
46 }
47 
48 
50 {
51 
52  // do anything here that needs to be done at destruction time
53  // (e.g. close files, deallocate resources etc.)
54 
55 }
56 
57 
58 //
59 // member functions
60 //
61 
62 // ------------ method called to for each event ------------
63 void
65 {
66  using namespace edm;
67 
68 #ifdef THIS_IS_AN_EVENT_EXAMPLE
70  iEvent.getByLabel("example",pIn);
71 #endif
72 
73 #ifdef THIS_IS_AN_EVENTSETUP_EXAMPLE
74  ESHandle<SetupData> pSetup;
75  iSetup.get<SetupRecord>().get(pSetup);
76 #endif
77 
78  ESHandle<L1CaloEcalScale> caloEcalScale;
79  ESHandle<L1CaloHcalScale> caloHcalScale;
80  ESHandle<CaloTPGTranscoder> caloTPGTranscoder;
81  iSetup.get<L1CaloEcalScaleRcd>().get(caloEcalScale);
82  iSetup.get<L1CaloHcalScaleRcd>().get(caloHcalScale);
83  iSetup.get<CaloTPGRecord>().get(caloTPGTranscoder);
84 
85  EcalTPGScale* ecalTPGScale = new EcalTPGScale();
86  ecalTPGScale->setEventSetup(iSetup);
87 
88  bool ecalIsConsistent = true;
89  bool hcalIsConsistent = true;
90 
91  double ecal1;
92  double ecal2;
93  double hcal1;
94  double hcal2;
95  double hcal3;
96  double hcal4;
97 
98 
99  // compare the ecal scales
100 
101  // 8 bits of input energy
102  for (unsigned short input = 0; input <= 0xFF; input++)
103  {
104  // loop over ietas, barrel first
105  for (unsigned short absIeta = 1; absIeta <= 17; absIeta++)
106  {
107 
108  // positive eta
109  ecal1 = ecalTPGScale->getTPGInGeV( (unsigned int) input,
111  absIeta, 1));
112  ecal2 = caloEcalScale->et(input, absIeta, 1);
113 
114  if ( !(ecal1 == ecal2) )
115  {
116  ecalIsConsistent = false;
117  /*LogWarning("InconsistentData")
118  << "ECAL scales not consistent! (pos eta barrel)"
119  << "old ECAL is " << setprecision (8) << ecal1
120  << " new ECAL is " << setprecision (8) << ecal2
121  << " difference is " << (ecal1 - ecal2) ; */
122  }
123 
124  // negative eta
125  ecal1 = ecalTPGScale->
126  getTPGInGeV( (unsigned int) input, EcalTrigTowerDetId(-1, EcalBarrel,
127  absIeta, 2));
128  ecal2 = caloEcalScale->et(input, absIeta, -1);
129 
130  if ( !(ecal1 == ecal2) )
131  {
132  ecalIsConsistent = false;
133  /*LogWarning("InconsistentData")
134  << "ECAL scales not consistent! (neg eta barrel)"
135  << "old ECAL is " << setprecision (8) << ecal1
136  << " new ECAL is " << setprecision (8) << ecal2
137  << " difference is " << (ecal1 - ecal2) ; */
138  }
139  }
140  // now loop over endcap ietas
141  for (unsigned short absIeta = 18; absIeta <= 28; absIeta++)
142  {
143 
144  // positive eta
145  ecal1 = ecalTPGScale->
146  getTPGInGeV( (unsigned int) input, EcalTrigTowerDetId(1, EcalEndcap,
147  absIeta, 1));
148  ecal2 = caloEcalScale->et(input, absIeta, 1);
149 
150  if ( !(ecal1 == ecal2) )
151  {
152  ecalIsConsistent = false;
153  /*LogWarning("InconsistentData")
154  << "ECAL scales not consistent! (pos eta endcap)"
155  << "old ECAL is " << setprecision (8) << ecal1
156  << " new ECAL is " << setprecision (8) << ecal2
157  << " difference is " << (ecal1 - ecal2) ; */
158  }
159 
160  // negative eta
161  ecal1 = ecalTPGScale->
162  getTPGInGeV( (unsigned int) input, EcalTrigTowerDetId(-1, EcalEndcap,
163  absIeta, 2));
164  ecal2 = caloEcalScale->et(input, absIeta, -1);
165 
166  if ( !(ecal1 == ecal2) )
167  {
168  ecalIsConsistent = false;
169  /*LogWarning("InconsistentData")
170  << "ECAL scales not consistent! (neg eta endcap)"
171  << "old ECAL is " << setprecision (8) << ecal1
172  << " new ECAL is " << setprecision (8) << ecal2
173  << " difference is " << (ecal1 - ecal2) ; */
174  }
175  }
176  }
177 
178  if (!ecalIsConsistent)
179  {
180  // do something
181  //cout << "WARNING: ECAL scales not consistent!" << endl;
182  LogWarning("InconsistentData") << "ECAL scales not consistent!";
183  }
184  else
185  {
186  // do something else
187  //cout << "ECAL scales okay" << endl;
188  }
189 
190  // compare the hcal scales
191 
192  for (unsigned short input = 0; input <= 0xFF; input++)
193  {
194  // loop over ietas
195  for (unsigned short absIeta = 1; absIeta <= 32; absIeta++)
196  {
197  hcal1 = caloTPGTranscoder->hcaletValue(absIeta, input); // no eta-
198  hcal2 = caloHcalScale->et(input, absIeta, 1); // sign in transcoder
199  hcal3 = caloTPGTranscoder->hcaletValue(-absIeta, input); // no eta-
200  hcal4 = caloHcalScale->et(input, absIeta,-1); // sign in transcoder
201 
202 
203 
204  if ( (!(hcal1 == hcal2))||(!(hcal3==hcal4)))
205  {
206  hcalIsConsistent = false;
207  /*LogWarning("InconsistentData")
208  << "HCAL scales not consistent!"
209  << "old HCAL is " << hcal1
210  << " new HCAL is " << hcal2 ; */
211  }
212  }
213  }
214  if (!hcalIsConsistent)
215  {
216  // do something
217  //cout << "WARNING: HCAL scales not consistent!" << endl;
218  LogWarning("InconsistentData") << "HCAL scales not consistent!";
219  }
220  else
221  {
222  // do something else
223  //cout << "HCAL scales okay" << endl;
224  }
225 }
226 
227 
228 // ------------ method called once each job just before starting event loop ------------
229 void
231 {
232 }
233 
234 // ------------ method called once each job just after ending the event loop ------------
235 void
237 }
238 
L1CaloInputScaleTester(const edm::ParameterSet &)
virtual void analyze(const edm::Event &, const edm::EventSetup &)
int iEvent
Definition: GenABIO.cc:243
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
const T & get() const
Definition: EventSetup.h:55
tuple cout
Definition: gather_cfg.py:121