CMS 3D CMS Logo

L1TPhysicalEtAdder.cc
Go to the documentation of this file.
3 
8 
12 
15 
18 
20 #include <vector>
21 
22 namespace {
23 
24 int getRegionEta(int gtEta, bool forward)
25 {
26  // backwards conversion is
27  // unsigned rctEta = (iEta<11 ? 10-iEta : iEta-11);
28  // return (((rctEta % 7) & 0x7) | (iEta<11 ? 0x8 : 0));
29  int centralGtEta[] = {11, 12, 13, 14, 15, 16, 17, -100, 10, 9, 8, 7, 6, 5, 4};
30  int forwardGtEta[] = {18, 19, 20, 21, -100, -100, -100, -100, 3, 2, 1, 0};
31 
32  //printf("%i, %i\n",gtEta,forward);
33 
34  int regionEta;
35 
36  if(!forward)
37  {
38  regionEta = centralGtEta[gtEta];
39  } else
40  regionEta = forwardGtEta[gtEta];
41 
42  if(regionEta == -100)
43  edm::LogError("EtaIndexError")
44  << "Bad eta index passed to L1TPhysicalEtAdder::getRegionEta, " << gtEta << std::endl;
45 
46  return regionEta;
47 }
48 
49 // adapted these from the UCT2015 codebase.
50 double getPhysicalEta(int gtEta, bool forward = false)
51 {
52  int etaIndex = getRegionEta(gtEta, forward);
53 
54  const double rgnEtaValues[11] = {
55  0.174, // HB and inner HE bins are 0.348 wide
56  0.522,
57  0.870,
58  1.218,
59  1.566,
60  1.956, // Last two HE bins are 0.432 and 0.828 wide
61  2.586,
62  3.250, // HF bins are 0.5 wide
63  3.750,
64  4.250,
65  4.750
66  };
67  if(etaIndex < 11) {
68  return -rgnEtaValues[-(etaIndex - 10)]; // 0-10 are negative eta values
69  }
70  else if (etaIndex < 22) {
71  return rgnEtaValues[etaIndex - 11]; // 11-21 are positive eta values
72  }
73  return -9;
74 }
75 
76 double getPhysicalPhi(int phiIndex)
77 {
78  if (phiIndex < 10)
79  return 2. * M_PI * phiIndex / 18.;
80  if (phiIndex < 18)
81  return -M_PI + 2. * M_PI * (phiIndex - 9) / 18.;
82  return -9;
83 }
84 
85 }
86 
87 using namespace l1t;
88 
90 
91  produces<EGammaBxCollection>();
92  produces<TauBxCollection>("rlxTaus");
93  produces<TauBxCollection>("isoTaus");
94  produces<JetBxCollection>();
95  produces<JetBxCollection>("preGtJets");
96  produces<EtSumBxCollection>();
97  produces<CaloSpareBxCollection>("HFRingSums");
98  produces<CaloSpareBxCollection>("HFBitCounts");
99 
100  EGammaToken_ = consumes<EGammaBxCollection>(ps.getParameter<edm::InputTag>("InputCollection"));
101  RlxTauToken_ = consumes<TauBxCollection>(ps.getParameter<edm::InputTag>("InputRlxTauCollection"));
102  IsoTauToken_ = consumes<TauBxCollection>(ps.getParameter<edm::InputTag>("InputIsoTauCollection"));
103  JetToken_ = consumes<JetBxCollection>(ps.getParameter<edm::InputTag>("InputCollection"));
104  preGtJetToken_ = consumes<JetBxCollection>(ps.getParameter<edm::InputTag>("InputPreGtJetCollection"));
105  EtSumToken_ = consumes<EtSumBxCollection>(ps.getParameter<edm::InputTag>("InputCollection"));
106  HfSumsToken_ = consumes<CaloSpareBxCollection>(ps.getParameter<edm::InputTag>("InputHFSumsCollection"));
107  HfCountsToken_ = consumes<CaloSpareBxCollection>(ps.getParameter<edm::InputTag>("InputHFCountsCollection"));
108 }
109 
111 
112 }
113 
114 // ------------ method called to produce the data ------------
115 void
117 {
118  // store new collections which include physical quantities
119  std::unique_ptr<EGammaBxCollection> new_egammas (new EGammaBxCollection);
120  std::unique_ptr<TauBxCollection> new_rlxtaus (new TauBxCollection);
121  std::unique_ptr<TauBxCollection> new_isotaus (new TauBxCollection);
122  std::unique_ptr<JetBxCollection> new_jets (new JetBxCollection);
123  std::unique_ptr<JetBxCollection> new_preGtJets (new JetBxCollection);
124  std::unique_ptr<EtSumBxCollection> new_etsums (new EtSumBxCollection);
125  std::unique_ptr<CaloSpareBxCollection> new_hfsums (new CaloSpareBxCollection);
126  std::unique_ptr<CaloSpareBxCollection> new_hfcounts (new CaloSpareBxCollection);
127 
129  edm::Handle<TauBxCollection> old_rlxtaus;
130  edm::Handle<TauBxCollection> old_isotaus;
132  edm::Handle<JetBxCollection> old_preGtJets;
136 
137  iEvent.getByToken(EGammaToken_, old_egammas);
138  iEvent.getByToken(RlxTauToken_, old_rlxtaus);
139  iEvent.getByToken(IsoTauToken_, old_isotaus);
140  iEvent.getByToken(JetToken_, old_jets);
141  iEvent.getByToken(preGtJetToken_, old_preGtJets);
142  iEvent.getByToken(EtSumToken_, old_etsums);
143  iEvent.getByToken(HfSumsToken_, old_hfsums);
144  iEvent.getByToken(HfCountsToken_, old_hfcounts);
145 
146  //get the proper scales for conversion to physical et
148  iSetup.get< L1EmEtScaleRcd >().get( emScale ) ;
149 
151  iSetup.get< L1JetEtScaleRcd >().get( jetScale ) ;
152 
153  edm::ESHandle< L1CaloEtScale > htMissScale ;
154  iSetup.get< L1HtMissScaleRcd >().get( htMissScale ) ;
155 
156  int firstBX = old_egammas->getFirstBX();
157  int lastBX = old_egammas->getLastBX();
158 
159  new_egammas->setBXRange(firstBX, lastBX);
160  new_rlxtaus->setBXRange(firstBX, lastBX);
161  new_isotaus->setBXRange(firstBX, lastBX);
162  new_jets->setBXRange(firstBX, lastBX);
163  new_preGtJets->setBXRange(firstBX, lastBX);
164  new_etsums->setBXRange(firstBX, lastBX);
165  new_hfsums->setBXRange(firstBX, lastBX);
166  new_hfcounts->setBXRange(firstBX, lastBX);
167 
168  for(int bx = firstBX; bx <= lastBX; ++bx)
169  {
170  for(EGammaBxCollection::const_iterator itEGamma = old_egammas->begin(bx);
171  itEGamma != old_egammas->end(bx); ++itEGamma)
172  {
173  //const double pt = itEGamma->hwPt() * emScale->linearLsb();
174  const double et = emScale->et( itEGamma->hwPt() );
175  const double eta = getPhysicalEta(itEGamma->hwEta());
176  const double phi = getPhysicalPhi(itEGamma->hwPhi());
177  math::PtEtaPhiMLorentzVector p4(et, eta, phi, 0);
178 
179  EGamma eg(*&p4, itEGamma->hwPt(),
180  itEGamma->hwEta(), itEGamma->hwPhi(),
181  itEGamma->hwQual(), itEGamma->hwIso());
182  new_egammas->push_back(bx, *&eg);
183 
184 
185  }
186 
187  for(TauBxCollection::const_iterator itTau = old_rlxtaus->begin(bx);
188  itTau != old_rlxtaus->end(bx); ++itTau)
189  {
190  // use the full-circle conversion to match l1extra, accounts for linearLsb and max value automatically
191  //const uint16_t rankPt = jetScale->rank((uint16_t)itTau->hwPt());
192  //const double et = jetScale->et( rankPt ) ;
193 
194  // or use the emScale to get finer-grained et
195  //const double et = itTau->hwPt() * emScale->linearLsb();
196 
197  // we are now already in the rankPt
198  const double et = jetScale->et( itTau->hwPt() );
199 
200  const double eta = getPhysicalEta(itTau->hwEta());
201  const double phi = getPhysicalPhi(itTau->hwPhi());
202  math::PtEtaPhiMLorentzVector p4(et, eta, phi, 0);
203 
204  Tau tau(*&p4, itTau->hwPt(),
205  itTau->hwEta(), itTau->hwPhi(),
206  itTau->hwQual(), itTau->hwIso());
207  new_rlxtaus->push_back(bx, *&tau);
208 
209  }
210 
211  for(TauBxCollection::const_iterator itTau = old_isotaus->begin(bx);
212  itTau != old_isotaus->end(bx); ++itTau)
213  {
214  // use the full-circle conversion to match l1extra, accounts for linearLsb and max value automatically
215  //const uint16_t rankPt = jetScale->rank((uint16_t)itTau->hwPt());
216  //const double et = jetScale->et( rankPt ) ;
217 
218  // or use the emScale to get finer-grained et
219  //const double et = itTau->hwPt() * emScale->linearLsb();
220 
221  // we are now already in the rankPt
222  const double et = jetScale->et( itTau->hwPt() );
223 
224  const double eta = getPhysicalEta(itTau->hwEta());
225  const double phi = getPhysicalPhi(itTau->hwPhi());
226  math::PtEtaPhiMLorentzVector p4(et, eta, phi, 0);
227 
228  Tau tau(*&p4, itTau->hwPt(),
229  itTau->hwEta(), itTau->hwPhi(),
230  itTau->hwQual(), itTau->hwIso());
231  new_isotaus->push_back(bx, *&tau);
232 
233  }
234 
235  for(JetBxCollection::const_iterator itJet = old_jets->begin(bx);
236  itJet != old_jets->end(bx); ++itJet)
237  {
238  // use the full-circle conversion to match l1extra, accounts for linearLsb and max value automatically
239  //const uint16_t rankPt = jetScale->rank((uint16_t)itJet->hwPt());
240  //const double et = jetScale->et( rankPt ) ;
241 
242  // or use the emScale to get finer-grained et
243  //const double et = itJet->hwPt() * emScale->linearLsb();
244 
245  // we are now already in the rankPt
246  const double et = jetScale->et( itJet->hwPt() );
247 
248  const bool forward = ((itJet->hwQual() & 0x2) != 0);
249  const double eta = getPhysicalEta(itJet->hwEta(), forward);
250  const double phi = getPhysicalPhi(itJet->hwPhi());
251  math::PtEtaPhiMLorentzVector p4(et, eta, phi, 0);
252 
253  Jet jet(*&p4, itJet->hwPt(),
254  itJet->hwEta(), itJet->hwPhi(),
255  itJet->hwQual());
256  new_jets->push_back(bx, *&jet);
257 
258  }
259 
260  for(JetBxCollection::const_iterator itJet = old_preGtJets->begin(bx);
261  itJet != old_preGtJets->end(bx); ++itJet)
262  {
263  // use the full-circle conversion to match l1extra, accounts for linearLsb and max value automatically
264  //const uint16_t rankPt = jetScale->rank((uint16_t)itJet->hwPt());
265  //const double et = jetScale->et( rankPt ) ;
266 
267  // or use the emScale to get finer-grained et
268  const double et = itJet->hwPt() * emScale->linearLsb();
269 
270  // we are now already in the rankPt
271  //const double et = jetScale->et( itJet->hwPt() );
272 
273  const bool forward = ((itJet->hwQual() & 0x2) != 0);
274  const double eta = getPhysicalEta(itJet->hwEta(), forward);
275  const double phi = getPhysicalPhi(itJet->hwPhi());
276  math::PtEtaPhiMLorentzVector p4(et, eta, phi, 0);
277 
278  Jet jet(*&p4, itJet->hwPt(),
279  itJet->hwEta(), itJet->hwPhi(),
280  itJet->hwQual());
281  new_preGtJets->push_back(bx, *&jet);
282 
283  }
284 
285 
286  for(EtSumBxCollection::const_iterator itEtSum = old_etsums->begin(bx);
287  itEtSum != old_etsums->end(bx); ++itEtSum)
288  {
289  double et = itEtSum->hwPt() * emScale->linearLsb();
290  //hack while we figure out the right scales
291  //double et = emScale->et( itEtSum->hwPt() );
292  const EtSum::EtSumType sumType = itEtSum->getType();
293 
294  const double eta = getPhysicalEta(itEtSum->hwEta());
295  double phi = getPhysicalPhi(itEtSum->hwPhi());
296  if(sumType == EtSum::EtSumType::kMissingHt){
297  et = htMissScale->et( itEtSum->hwPt() );
298  double regionPhiWidth=2. * 3.1415927 / L1CaloRegionDetId::N_PHI;
299  phi=phi+(regionPhiWidth/2.); // add the region half-width to match L1Extra MHT phi
300  }
301 
302 
303  math::PtEtaPhiMLorentzVector p4(et, eta, phi, 0);
304 
305  EtSum eg(*&p4, sumType, itEtSum->hwPt(),
306  itEtSum->hwEta(), itEtSum->hwPhi(),
307  itEtSum->hwQual());
308  new_etsums->push_back(bx, *&eg);
309 
310 
311  }
312 
313  for(CaloSpareBxCollection::const_iterator itCaloSpare = old_hfsums->begin(bx);
314  itCaloSpare != old_hfsums->end(bx); ++itCaloSpare)
315  {
316  //just pass through for now
317  //a different scale is needed depending on the type
318  new_hfsums->push_back(bx, *itCaloSpare);
319  }
320 
321  for(CaloSpareBxCollection::const_iterator itCaloSpare = old_hfcounts->begin(bx);
322  itCaloSpare != old_hfcounts->end(bx); ++itCaloSpare)
323  {
324  //just pass through for now
325  //a different scale is needed depending on the type
326  new_hfcounts->push_back(bx, *itCaloSpare);
327  }
328 
329  }
330 
331  iEvent.put(std::move(new_egammas));
332  iEvent.put(std::move(new_rlxtaus),"rlxTaus");
333  iEvent.put(std::move(new_isotaus),"isoTaus");
334  iEvent.put(std::move(new_jets));
335  iEvent.put(std::move(new_preGtJets),"preGtJets");
336  iEvent.put(std::move(new_etsums));
337  iEvent.put(std::move(new_hfsums),"HFRingSums");
338  iEvent.put(std::move(new_hfcounts),"HFBitCounts");
339 }
340 
341 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
342 void
344  //The following says we do not know what parameters are allowed so do no validation
345  // Please change this to state exactly what you do use, even if it is no parameters
347  desc.setUnknown();
348  descriptions.addDefault(desc);
349 }
350 
351 //define this as a plug-in
353 
const_iterator end(int bx) const
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const unsigned int gtEta(const unsigned int iEta)
delete x;
Definition: CaloConfig.h:22
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
double et(const uint16_t rank) const
convert from rank to physically meaningful quantity
double linearLsb() const
get LSB of linear input scale
Definition: L1CaloEtScale.h:53
L1TPhysicalEtAdder(const edm::ParameterSet &ps)
PtEtaPhiMLorentzVectorD PtEtaPhiMLorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:25
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
Definition: Jet.py:1
double p4[4]
Definition: TauolaWrapper.h:92
#define M_PI
Definition: Tau.py:1
int getFirstBX() const
const T & get() const
Definition: EventSetup.h:58
et
define resolution functions of each parameter
int getLastBX() const
static const unsigned N_PHI
const_iterator begin(int bx) const
EtSumType
Definition: EtSum.h:20
def move(src, dest)
Definition: eostools.py:510
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< EGamma >::const_iterator const_iterator
Definition: BXVector.h:20