CMS 3D CMS Logo

L1CaloInputScalesGenerator.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 #include <fstream>
11 using std::ofstream;
12 
13 // user include files
15 
19 
21 
26 
27 //
28 // constants, enums and typedefs
29 //
30 
31 //
32 // static data member definitions
33 //
34 
35 //
36 // constructors and destructor
37 //
39 
40 {
41  //now do what ever initialization is needed
42 }
43 
45  // do anything here that needs to be done at desctruction time
46  // (e.g. close files, deallocate resources etc.)
47 }
48 
49 //
50 // member functions
51 //
52 
53 // ------------ method called to for each event ------------
55  using namespace edm;
56 
57  ESHandle<CaloTPGTranscoder> caloTPGTranscoder;
58  iSetup.get<CaloTPGRecord>().get(caloTPGTranscoder);
59 
60  EcalTPGScale* ecalTPGScale = new EcalTPGScale();
61  ecalTPGScale->setEventSetup(iSetup);
62 
63  double output;
64  ofstream scalesFile("L1CaloInputScales_cfi.py");
65 
66  // Write the ecal scales, positive eta
67 
68  scalesFile << "import FWCore.ParameterSet.Config as cms\n" << endl;
69 
70  scalesFile << "L1CaloInputScalesProducer =cms.ESProducer(\"L1CaloInputScalesProducer\"," << endl;
71  scalesFile << "L1EcalEtThresholdsPositiveEta = cms.vdouble(" << endl;
72 
73  //Python does not support arrays over 255 entries so we neeed ton accomodate it by creating new array after 255 entries
74  int nEntries = 0;
75 
76  // loop over ietas, barrel
77  for (unsigned short absIeta = 1; absIeta <= 28; absIeta++) {
78  EcalSubdetector subdet = (absIeta <= 17) ? EcalBarrel : EcalEndcap;
79  // 8 bits of input energy
80  for (unsigned short input = 0; input <= 0xFF; input++) {
81  output = ecalTPGScale->getTPGInGeV((unsigned int)input, EcalTrigTowerDetId(1, subdet, absIeta, 1));
82  scalesFile << setprecision(8) << output;
83  nEntries++;
84 
85  if (absIeta == 28 && input == 0xFF) {
86  scalesFile << "),";
87  } else if (nEntries > 254) {
88  scalesFile << ")+cms.vdouble(";
89  nEntries = 0;
90  } else {
91  scalesFile << ", ";
92  }
93  }
94  scalesFile << endl;
95  }
96 
97  // Write the ecal scales, negative eta
98 
99  scalesFile << endl << "\tL1EcalEtThresholdsNegativeEta = cms.vdouble(" << endl;
100 
101  nEntries = 0;
102  // loop over ietas, barrel first
103  for (unsigned short absIeta = 1; absIeta <= 28; absIeta++) {
104  EcalSubdetector subdet = (absIeta <= 17) ? EcalBarrel : EcalEndcap;
105  // 8 bits of input energy
106  for (unsigned short input = 0; input <= 0xFF; input++) {
107  // negative eta
108  output = ecalTPGScale->getTPGInGeV((unsigned int)input, EcalTrigTowerDetId(-1, subdet, absIeta, 2));
109  scalesFile << setprecision(8) << output;
110  nEntries++;
111 
112  if (absIeta == 28 && input == 0xFF) {
113  scalesFile << "),";
114  } else if (nEntries > 254) {
115  scalesFile << ")+cms.vdouble(";
116  nEntries = 0;
117  } else {
118  scalesFile << ", ";
119  }
120  }
121  scalesFile << endl;
122  }
123 
124  // Write the hcal scales (Positive Eta)
125 
126  scalesFile << endl << "\tL1HcalEtThresholdsPositiveEta = cms.vdouble(" << endl;
127 
128  // loop over ietas
129 
130  nEntries = 0;
131  for (unsigned short absIeta = 1; absIeta <= 32; absIeta++) {
132  for (unsigned short input = 0; input <= 0xFF; input++) {
133  output = caloTPGTranscoder->hcaletValue(absIeta, input);
134  scalesFile << setprecision(8) << output;
135  nEntries++;
136 
137  if (absIeta == 32 && input == 0xFF) {
138  scalesFile << "),";
139  } else if (nEntries > 254) {
140  scalesFile << ")+cms.vdouble(";
141  nEntries = 0;
142  } else {
143  scalesFile << ", ";
144  }
145  }
146  scalesFile << endl;
147  }
148 
149  // Write the hcal scales (Negative Eta)
150 
151  scalesFile << endl << "\tL1HcalEtThresholdsNegativeEta = cms.vdouble(" << endl;
152 
153  nEntries = 0;
154  // loop over ietas
155  for (unsigned short absIeta = 1; absIeta <= 32; absIeta++) {
156  for (unsigned short input = 0; input <= 0xFF; input++) {
157  output = caloTPGTranscoder->hcaletValue(-absIeta, input);
158  scalesFile << setprecision(8) << output;
159  nEntries++;
160 
161  if (absIeta == 32 && input == 0xFF) {
162  scalesFile << ")";
163  } else if (nEntries > 254) {
164  scalesFile << ")+cms.vdouble(";
165  nEntries = 0;
166  } else {
167  scalesFile << ", ";
168  }
169  }
170  scalesFile << endl;
171  }
172 
173  scalesFile << ")" << endl;
174 
175  scalesFile.close();
176 }
177 
178 // ------------ method called once each job just before starting event loop ------------
180 
181 // ------------ method called once each job just after ending the event loop ------------
EcalTPGScale::getTPGInGeV
double getTPGInGeV(const EcalTriggerPrimitiveDigi &tpDigi)
Definition: EcalTPGScale.cc:18
input
static const std::string input
Definition: EdmProvDump.cc:48
EcalTPGScale
Definition: EcalTPGScale.h:8
MessageLogger.h
ESHandle.h
convertSQLitetoXML_cfg.output
output
Definition: convertSQLitetoXML_cfg.py:32
edm
HLT enums.
Definition: AlignableModifier.h:19
gather_cfg.cout
cout
Definition: gather_cfg.py:144
L1CaloInputScalesGenerator::endJob
void endJob() override
Definition: L1CaloInputScalesGenerator.cc:182
EcalSubdetector
EcalSubdetector
Definition: EcalSubdetector.h:10
EcalTPGScale::setEventSetup
void setEventSetup(const edm::EventSetup &evtSetup)
Definition: EcalTPGScale.cc:16
EcalTrigTowerDetId
Definition: EcalTrigTowerDetId.h:14
CaloTPGTranscoder.h
EcalBarrel
Definition: EcalSubdetector.h:10
MakerMacros.h
edm::EventSetup::get
T get() const
Definition: EventSetup.h:73
edm::ESHandle< CaloTPGTranscoder >
EcalTrigTowerDetId.h
EcalEndcap
Definition: EcalSubdetector.h:10
L1CaloInputScalesGenerator.h
L1CaloInputScalesGenerator::~L1CaloInputScalesGenerator
~L1CaloInputScalesGenerator() override
Definition: L1CaloInputScalesGenerator.cc:44
edm::ParameterSet
Definition: ParameterSet.h:36
CaloTPGRecord
Definition: CaloTPGRecord.h:26
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::EventSetup
Definition: EventSetup.h:57
CaloTPGRecord.h
get
#define get
L1CaloInputScalesGenerator::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: L1CaloInputScalesGenerator.cc:54
Frameworkfwd.h
EventSetup.h
CaloTPGTranscoder::hcaletValue
virtual double hcaletValue(const int &ieta, const int &iphi, const int &version, const int &compressedValue) const =0
L1CaloInputScalesGenerator::L1CaloInputScalesGenerator
L1CaloInputScalesGenerator(const edm::ParameterSet &)
Definition: L1CaloInputScalesGenerator.cc:38
L1CaloInputScalesGenerator::beginJob
void beginJob() override
Definition: L1CaloInputScalesGenerator.cc:179
edm::Event
Definition: Event.h:73
EcalTPGScale.h