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