Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-10-29 06:08:30

0001 #include "L1Trigger/L1TCalorimeter/plugins/L1TPhysicalEtAdder.h"
0002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0003 
0004 #include "DataFormats/L1CaloTrigger/interface/L1CaloEmCand.h"
0005 #include "DataFormats/L1CaloTrigger/interface/L1CaloRegion.h"
0006 #include "DataFormats/L1CaloTrigger/interface/L1CaloCollections.h"
0007 
0008 #include "DataFormats/L1Trigger/interface/BXVector.h"
0009 #include "DataFormats/L1Trigger/interface/EtSum.h"
0010 
0011 #include "DataFormats/L1TCalorimeter/interface/CaloEmCand.h"
0012 #include "DataFormats/L1TCalorimeter/interface/CaloRegion.h"
0013 
0014 #include "FWCore/Framework/interface/ESHandle.h"
0015 #include <cassert>
0016 #include <vector>
0017 
0018 namespace {
0019 
0020   int getRegionEta(int gtEta, bool forward) {
0021     // backwards conversion is
0022     // unsigned rctEta = (iEta<11 ? 10-iEta : iEta-11);
0023     // return (((rctEta % 7) & 0x7) | (iEta<11 ? 0x8 : 0));
0024     int centralGtEta[] = {11, 12, 13, 14, 15, 16, 17, -100, 10, 9, 8, 7, 6, 5, 4};
0025     int forwardGtEta[] = {18, 19, 20, 21, -100, -100, -100, -100, 3, 2, 1, 0};
0026 
0027     //printf("%i, %i\n",gtEta,forward);
0028 
0029     int regionEta;
0030 
0031     if (!forward) {
0032       regionEta = centralGtEta[gtEta];
0033     } else
0034       regionEta = forwardGtEta[gtEta];
0035 
0036     if (regionEta == -100)
0037       edm::LogError("EtaIndexError") << "Bad eta index passed to L1TPhysicalEtAdder::getRegionEta, " << gtEta
0038                                      << std::endl;
0039 
0040     return regionEta;
0041   }
0042 
0043   // adapted these from the UCT2015 codebase.
0044   double getPhysicalEta(int gtEta, bool forward = false) {
0045     int etaIndex = getRegionEta(gtEta, forward);
0046 
0047     const double rgnEtaValues[11] = {0.174,  // HB and inner HE bins are 0.348 wide
0048                                      0.522,
0049                                      0.870,
0050                                      1.218,
0051                                      1.566,
0052                                      1.956,  // Last two HE bins are 0.432 and 0.828 wide
0053                                      2.586,
0054                                      3.250,  // HF bins are 0.5 wide
0055                                      3.750,
0056                                      4.250,
0057                                      4.750};
0058     if (0 <= etaIndex && etaIndex < 11) {
0059       return -rgnEtaValues[-(etaIndex - 10)];  // 0-10 are negative eta values
0060     } else if (etaIndex < 22) {
0061       return rgnEtaValues[etaIndex - 11];  // 11-21 are positive eta values
0062     }
0063     return -9;
0064   }
0065 
0066   double getPhysicalPhi(int phiIndex) {
0067     if (phiIndex < 10)
0068       return 2. * M_PI * phiIndex / 18.;
0069     if (phiIndex < 18)
0070       return -M_PI + 2. * M_PI * (phiIndex - 9) / 18.;
0071     return -9;
0072   }
0073 
0074 }  // namespace
0075 
0076 using namespace l1t;
0077 
0078 L1TPhysicalEtAdder::L1TPhysicalEtAdder(const edm::ParameterSet& ps) {
0079   produces<EGammaBxCollection>();
0080   produces<TauBxCollection>("rlxTaus");
0081   produces<TauBxCollection>("isoTaus");
0082   produces<JetBxCollection>();
0083   produces<JetBxCollection>("preGtJets");
0084   produces<EtSumBxCollection>();
0085   produces<CaloSpareBxCollection>("HFRingSums");
0086   produces<CaloSpareBxCollection>("HFBitCounts");
0087 
0088   EGammaToken_ = consumes<EGammaBxCollection>(ps.getParameter<edm::InputTag>("InputCollection"));
0089   RlxTauToken_ = consumes<TauBxCollection>(ps.getParameter<edm::InputTag>("InputRlxTauCollection"));
0090   IsoTauToken_ = consumes<TauBxCollection>(ps.getParameter<edm::InputTag>("InputIsoTauCollection"));
0091   JetToken_ = consumes<JetBxCollection>(ps.getParameter<edm::InputTag>("InputCollection"));
0092   preGtJetToken_ = consumes<JetBxCollection>(ps.getParameter<edm::InputTag>("InputPreGtJetCollection"));
0093   EtSumToken_ = consumes<EtSumBxCollection>(ps.getParameter<edm::InputTag>("InputCollection"));
0094   HfSumsToken_ = consumes<CaloSpareBxCollection>(ps.getParameter<edm::InputTag>("InputHFSumsCollection"));
0095   HfCountsToken_ = consumes<CaloSpareBxCollection>(ps.getParameter<edm::InputTag>("InputHFCountsCollection"));
0096   emScaleToken_ = esConsumes<L1CaloEtScale, L1EmEtScaleRcd>();
0097   jetScaleToken_ = esConsumes<L1CaloEtScale, L1JetEtScaleRcd>();
0098   htMissScaleToken_ = esConsumes<L1CaloEtScale, L1HtMissScaleRcd>();
0099 }
0100 
0101 L1TPhysicalEtAdder::~L1TPhysicalEtAdder() {}
0102 
0103 // ------------ method called to produce the data  ------------
0104 void L1TPhysicalEtAdder::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0105   // store new collections which include physical quantities
0106   std::unique_ptr<EGammaBxCollection> new_egammas(new EGammaBxCollection);
0107   std::unique_ptr<TauBxCollection> new_rlxtaus(new TauBxCollection);
0108   std::unique_ptr<TauBxCollection> new_isotaus(new TauBxCollection);
0109   std::unique_ptr<JetBxCollection> new_jets(new JetBxCollection);
0110   std::unique_ptr<JetBxCollection> new_preGtJets(new JetBxCollection);
0111   std::unique_ptr<EtSumBxCollection> new_etsums(new EtSumBxCollection);
0112   std::unique_ptr<CaloSpareBxCollection> new_hfsums(new CaloSpareBxCollection);
0113   std::unique_ptr<CaloSpareBxCollection> new_hfcounts(new CaloSpareBxCollection);
0114 
0115   edm::Handle<EGammaBxCollection> old_egammas;
0116   edm::Handle<TauBxCollection> old_rlxtaus;
0117   edm::Handle<TauBxCollection> old_isotaus;
0118   edm::Handle<JetBxCollection> old_jets;
0119   edm::Handle<JetBxCollection> old_preGtJets;
0120   edm::Handle<EtSumBxCollection> old_etsums;
0121   edm::Handle<CaloSpareBxCollection> old_hfsums;
0122   edm::Handle<CaloSpareBxCollection> old_hfcounts;
0123 
0124   iEvent.getByToken(EGammaToken_, old_egammas);
0125   iEvent.getByToken(RlxTauToken_, old_rlxtaus);
0126   iEvent.getByToken(IsoTauToken_, old_isotaus);
0127   iEvent.getByToken(JetToken_, old_jets);
0128   iEvent.getByToken(preGtJetToken_, old_preGtJets);
0129   iEvent.getByToken(EtSumToken_, old_etsums);
0130   iEvent.getByToken(HfSumsToken_, old_hfsums);
0131   iEvent.getByToken(HfCountsToken_, old_hfcounts);
0132 
0133   //get the proper scales for conversion to physical et
0134   edm::ESHandle<L1CaloEtScale> emScale = iSetup.getHandle(emScaleToken_);
0135 
0136   edm::ESHandle<L1CaloEtScale> jetScale = iSetup.getHandle(jetScaleToken_);
0137 
0138   edm::ESHandle<L1CaloEtScale> htMissScale = iSetup.getHandle(htMissScaleToken_);
0139 
0140   int firstBX = old_egammas->getFirstBX();
0141   int lastBX = old_egammas->getLastBX();
0142 
0143   new_egammas->setBXRange(firstBX, lastBX);
0144   new_rlxtaus->setBXRange(firstBX, lastBX);
0145   new_isotaus->setBXRange(firstBX, lastBX);
0146   new_jets->setBXRange(firstBX, lastBX);
0147   new_preGtJets->setBXRange(firstBX, lastBX);
0148   new_etsums->setBXRange(firstBX, lastBX);
0149   new_hfsums->setBXRange(firstBX, lastBX);
0150   new_hfcounts->setBXRange(firstBX, lastBX);
0151 
0152   for (int bx = firstBX; bx <= lastBX; ++bx) {
0153     for (EGammaBxCollection::const_iterator itEGamma = old_egammas->begin(bx); itEGamma != old_egammas->end(bx);
0154          ++itEGamma) {
0155       //const double pt = itEGamma->hwPt() * emScale->linearLsb();
0156       const double et = emScale->et(itEGamma->hwPt());
0157       const double eta = getPhysicalEta(itEGamma->hwEta());
0158       const double phi = getPhysicalPhi(itEGamma->hwPhi());
0159       math::PtEtaPhiMLorentzVector p4(et, eta, phi, 0);
0160 
0161       EGamma eg(*&p4, itEGamma->hwPt(), itEGamma->hwEta(), itEGamma->hwPhi(), itEGamma->hwQual(), itEGamma->hwIso());
0162       new_egammas->push_back(bx, *&eg);
0163     }
0164 
0165     for (TauBxCollection::const_iterator itTau = old_rlxtaus->begin(bx); itTau != old_rlxtaus->end(bx); ++itTau) {
0166       // use the full-circle conversion to match l1extra, accounts for linearLsb and max value automatically
0167       //const uint16_t rankPt = jetScale->rank((uint16_t)itTau->hwPt());
0168       //const double et = jetScale->et( rankPt ) ;
0169 
0170       // or use the emScale to get finer-grained et
0171       //const double et = itTau->hwPt() * emScale->linearLsb();
0172 
0173       // we are now already in the rankPt
0174       const double et = jetScale->et(itTau->hwPt());
0175 
0176       const double eta = getPhysicalEta(itTau->hwEta());
0177       const double phi = getPhysicalPhi(itTau->hwPhi());
0178       math::PtEtaPhiMLorentzVector p4(et, eta, phi, 0);
0179 
0180       Tau tau(*&p4, itTau->hwPt(), itTau->hwEta(), itTau->hwPhi(), itTau->hwQual(), itTau->hwIso());
0181       new_rlxtaus->push_back(bx, *&tau);
0182     }
0183 
0184     for (TauBxCollection::const_iterator itTau = old_isotaus->begin(bx); itTau != old_isotaus->end(bx); ++itTau) {
0185       // use the full-circle conversion to match l1extra, accounts for linearLsb and max value automatically
0186       //const uint16_t rankPt = jetScale->rank((uint16_t)itTau->hwPt());
0187       //const double et = jetScale->et( rankPt ) ;
0188 
0189       // or use the emScale to get finer-grained et
0190       //const double et = itTau->hwPt() * emScale->linearLsb();
0191 
0192       // we are now already in the rankPt
0193       const double et = jetScale->et(itTau->hwPt());
0194 
0195       const double eta = getPhysicalEta(itTau->hwEta());
0196       const double phi = getPhysicalPhi(itTau->hwPhi());
0197       math::PtEtaPhiMLorentzVector p4(et, eta, phi, 0);
0198 
0199       Tau tau(*&p4, itTau->hwPt(), itTau->hwEta(), itTau->hwPhi(), itTau->hwQual(), itTau->hwIso());
0200       new_isotaus->push_back(bx, *&tau);
0201     }
0202 
0203     for (JetBxCollection::const_iterator itJet = old_jets->begin(bx); itJet != old_jets->end(bx); ++itJet) {
0204       // use the full-circle conversion to match l1extra, accounts for linearLsb and max value automatically
0205       //const uint16_t rankPt = jetScale->rank((uint16_t)itJet->hwPt());
0206       //const double et = jetScale->et( rankPt ) ;
0207 
0208       // or use the emScale to get finer-grained et
0209       //const double et = itJet->hwPt() * emScale->linearLsb();
0210 
0211       // we are now already in the rankPt
0212       const double et = jetScale->et(itJet->hwPt());
0213 
0214       const bool forward = ((itJet->hwQual() & 0x2) != 0);
0215       const double eta = getPhysicalEta(itJet->hwEta(), forward);
0216       const double phi = getPhysicalPhi(itJet->hwPhi());
0217       math::PtEtaPhiMLorentzVector p4(et, eta, phi, 0);
0218 
0219       Jet jet(*&p4, itJet->hwPt(), itJet->hwEta(), itJet->hwPhi(), itJet->hwQual());
0220       new_jets->push_back(bx, *&jet);
0221     }
0222 
0223     for (JetBxCollection::const_iterator itJet = old_preGtJets->begin(bx); itJet != old_preGtJets->end(bx); ++itJet) {
0224       // use the full-circle conversion to match l1extra, accounts for linearLsb and max value automatically
0225       //const uint16_t rankPt = jetScale->rank((uint16_t)itJet->hwPt());
0226       //const double et = jetScale->et( rankPt ) ;
0227 
0228       // or use the emScale to get finer-grained et
0229       const double et = itJet->hwPt() * emScale->linearLsb();
0230 
0231       // we are now already in the rankPt
0232       //const double et = jetScale->et( itJet->hwPt() );
0233 
0234       const bool forward = ((itJet->hwQual() & 0x2) != 0);
0235       const double eta = getPhysicalEta(itJet->hwEta(), forward);
0236       const double phi = getPhysicalPhi(itJet->hwPhi());
0237       math::PtEtaPhiMLorentzVector p4(et, eta, phi, 0);
0238 
0239       Jet jet(*&p4, itJet->hwPt(), itJet->hwEta(), itJet->hwPhi(), itJet->hwQual());
0240       new_preGtJets->push_back(bx, *&jet);
0241     }
0242 
0243     for (EtSumBxCollection::const_iterator itEtSum = old_etsums->begin(bx); itEtSum != old_etsums->end(bx); ++itEtSum) {
0244       double et = itEtSum->hwPt() * emScale->linearLsb();
0245       //hack while we figure out the right scales
0246       //double et = emScale->et( itEtSum->hwPt() );
0247       const EtSum::EtSumType sumType = itEtSum->getType();
0248 
0249       const double eta = getPhysicalEta(itEtSum->hwEta());
0250       double phi = getPhysicalPhi(itEtSum->hwPhi());
0251       if (sumType == EtSum::EtSumType::kMissingHt) {
0252         et = htMissScale->et(itEtSum->hwPt());
0253         double regionPhiWidth = 2. * 3.1415927 / L1CaloRegionDetId::N_PHI;
0254         phi = phi + (regionPhiWidth / 2.);  // add the region half-width to match L1Extra MHT phi
0255       }
0256 
0257       math::PtEtaPhiMLorentzVector p4(et, eta, phi, 0);
0258 
0259       EtSum eg(*&p4, sumType, itEtSum->hwPt(), itEtSum->hwEta(), itEtSum->hwPhi(), itEtSum->hwQual());
0260       new_etsums->push_back(bx, *&eg);
0261     }
0262 
0263     for (CaloSpareBxCollection::const_iterator itCaloSpare = old_hfsums->begin(bx); itCaloSpare != old_hfsums->end(bx);
0264          ++itCaloSpare) {
0265       //just pass through for now
0266       //a different scale is needed depending on the type
0267       new_hfsums->push_back(bx, *itCaloSpare);
0268     }
0269 
0270     for (CaloSpareBxCollection::const_iterator itCaloSpare = old_hfcounts->begin(bx);
0271          itCaloSpare != old_hfcounts->end(bx);
0272          ++itCaloSpare) {
0273       //just pass through for now
0274       //a different scale is needed depending on the type
0275       new_hfcounts->push_back(bx, *itCaloSpare);
0276     }
0277   }
0278 
0279   iEvent.put(std::move(new_egammas));
0280   iEvent.put(std::move(new_rlxtaus), "rlxTaus");
0281   iEvent.put(std::move(new_isotaus), "isoTaus");
0282   iEvent.put(std::move(new_jets));
0283   iEvent.put(std::move(new_preGtJets), "preGtJets");
0284   iEvent.put(std::move(new_etsums));
0285   iEvent.put(std::move(new_hfsums), "HFRingSums");
0286   iEvent.put(std::move(new_hfcounts), "HFBitCounts");
0287 }
0288 
0289 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0290 void L1TPhysicalEtAdder::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0291   //The following says we do not know what parameters are allowed so do no validation
0292   // Please change this to state exactly what you do use, even if it is no parameters
0293   edm::ParameterSetDescription desc;
0294   desc.setUnknown();
0295   descriptions.addDefault(desc);
0296 }
0297 
0298 //define this as a plug-in
0299 DEFINE_FWK_MODULE(L1TPhysicalEtAdder);