Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:22:22

0001 #include "L1TriggerConfig/L1ScalesProducers/interface/L1CaloInputScaleTester.h"
0002 
0003 // system include files
0004 #include <memory>
0005 #include <iostream>
0006 using std::cout;
0007 using std::endl;
0008 #include <iomanip>
0009 using std::setprecision;
0010 
0011 // user include files
0012 #include "FWCore/Framework/interface/Frameworkfwd.h"
0013 
0014 #include "FWCore/Framework/interface/MakerMacros.h"
0015 #include "FWCore/Framework/interface/EventSetup.h"
0016 #include "FWCore/Framework/interface/ESHandle.h"
0017 
0018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0019 
0020 #include "CalibFormats/CaloTPG/interface/CaloTPGTranscoder.h"
0021 #include "CalibFormats/CaloTPG/interface/CaloTPGRecord.h"
0022 #include "CalibCalorimetry/EcalTPGTools/interface/EcalTPGScale.h"
0023 #include "CondFormats/L1TObjects/interface/L1CaloEcalScale.h"
0024 #include "CondFormats/DataRecord/interface/L1CaloEcalScaleRcd.h"
0025 #include "CondFormats/L1TObjects/interface/L1CaloHcalScale.h"
0026 #include "CondFormats/DataRecord/interface/L1CaloHcalScaleRcd.h"
0027 #include "DataFormats/EcalDetId/interface/EcalTrigTowerDetId.h"
0028 
0029 //
0030 // constants, enums and typedefs
0031 //
0032 
0033 //
0034 // static data member definitions
0035 //
0036 
0037 //
0038 // constructors and destructor
0039 //
0040 L1CaloInputScaleTester::L1CaloInputScaleTester(const edm::ParameterSet& iConfig)
0041     : hcalScaleToken_(esConsumes<L1CaloHcalScale, L1CaloHcalScaleRcd>()),
0042       ecalScaleToken_(esConsumes<L1CaloEcalScale, L1CaloEcalScaleRcd>()),
0043       transcoderToken_(esConsumes<CaloTPGTranscoder, CaloTPGRecord>()),
0044       tokens_(consumesCollector()) {
0045   //now do what ever initialization is needed
0046 }
0047 
0048 L1CaloInputScaleTester::~L1CaloInputScaleTester() {
0049   // do anything here that needs to be done at destruction time
0050   // (e.g. close files, deallocate resources etc.)
0051 }
0052 
0053 //
0054 // member functions
0055 //
0056 
0057 // ------------ method called to for each event  ------------
0058 void L1CaloInputScaleTester::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0059   using namespace edm;
0060 
0061 #ifdef THIS_IS_AN_EVENT_EXAMPLE
0062   Handle<ExampleData> pIn;
0063   iEvent.getByLabel("example", pIn);
0064 #endif
0065 
0066 #ifdef THIS_IS_AN_EVENTSETUP_EXAMPLE
0067   ESHandle<SetupData> pSetup;
0068   iSetup.get<SetupRecord>().get(pSetup);
0069 #endif
0070 
0071   ESHandle<L1CaloEcalScale> caloEcalScale = iSetup.getHandle(ecalScaleToken_);
0072   ESHandle<L1CaloHcalScale> caloHcalScale = iSetup.getHandle(hcalScaleToken_);
0073   ESHandle<CaloTPGTranscoder> caloTPGTranscoder = iSetup.getHandle(transcoderToken_);
0074 
0075   EcalTPGScale ecalTPGScale(tokens_, iSetup);
0076 
0077   bool ecalIsConsistent = true;
0078   bool hcalIsConsistent = true;
0079 
0080   double ecal1;
0081   double ecal2;
0082   double hcal1;
0083   double hcal2;
0084   double hcal3;
0085   double hcal4;
0086 
0087   // compare the ecal scales
0088 
0089   // 8 bits of input energy
0090   for (unsigned short input = 0; input <= 0xFF; input++) {
0091     // loop over ietas, barrel first
0092     for (unsigned short absIeta = 1; absIeta <= 17; absIeta++) {
0093       // positive eta
0094       ecal1 = ecalTPGScale.getTPGInGeV((unsigned int)input, EcalTrigTowerDetId(1, EcalBarrel, absIeta, 1));
0095       ecal2 = caloEcalScale->et(input, absIeta, 1);
0096 
0097       if (!(ecal1 == ecal2)) {
0098         ecalIsConsistent = false;
0099         /*LogWarning("InconsistentData") 
0100          << "ECAL scales not consistent! (pos eta barrel)"
0101          << "old ECAL is " << setprecision (8) << ecal1
0102          << " new ECAL is " << setprecision (8) << ecal2 
0103          << " difference is " << (ecal1 - ecal2) ; */
0104       }
0105 
0106       // negative eta
0107       ecal1 = ecalTPGScale.getTPGInGeV((unsigned int)input, EcalTrigTowerDetId(-1, EcalBarrel, absIeta, 2));
0108       ecal2 = caloEcalScale->et(input, absIeta, -1);
0109 
0110       if (!(ecal1 == ecal2)) {
0111         ecalIsConsistent = false;
0112         /*LogWarning("InconsistentData") 
0113          << "ECAL scales not consistent! (neg eta barrel)"
0114          << "old ECAL is " << setprecision (8) << ecal1
0115          << " new ECAL is " << setprecision (8) << ecal2           
0116          << " difference is " << (ecal1 - ecal2) ; */
0117       }
0118     }
0119     // now loop over endcap ietas
0120     for (unsigned short absIeta = 18; absIeta <= 28; absIeta++) {
0121       // positive eta
0122       ecal1 = ecalTPGScale.getTPGInGeV((unsigned int)input, EcalTrigTowerDetId(1, EcalEndcap, absIeta, 1));
0123       ecal2 = caloEcalScale->et(input, absIeta, 1);
0124 
0125       if (!(ecal1 == ecal2)) {
0126         ecalIsConsistent = false;
0127         /*LogWarning("InconsistentData") 
0128          << "ECAL scales not consistent! (pos eta endcap)"
0129          << "old ECAL is " << setprecision (8) << ecal1
0130          << " new ECAL is " << setprecision (8) << ecal2 
0131          << " difference is " << (ecal1 - ecal2) ; */
0132       }
0133 
0134       // negative eta
0135       ecal1 = ecalTPGScale.getTPGInGeV((unsigned int)input, EcalTrigTowerDetId(-1, EcalEndcap, absIeta, 2));
0136       ecal2 = caloEcalScale->et(input, absIeta, -1);
0137 
0138       if (!(ecal1 == ecal2)) {
0139         ecalIsConsistent = false;
0140         /*LogWarning("InconsistentData") 
0141          << "ECAL scales not consistent! (neg eta endcap)"
0142          << "old ECAL is " << setprecision (8) << ecal1
0143          << " new ECAL is " << setprecision (8) << ecal2 
0144          << " difference is " << (ecal1 - ecal2) ; */
0145       }
0146     }
0147   }
0148 
0149   if (!ecalIsConsistent) {
0150     // do something
0151     //cout << "WARNING: ECAL scales not consistent!" << endl;
0152     LogWarning("InconsistentData") << "ECAL scales not consistent!";
0153   } else {
0154     // do something else
0155     //cout << "ECAL scales okay" << endl;
0156   }
0157 
0158   // compare the hcal scales
0159 
0160   for (unsigned short input = 0; input <= 0xFF; input++) {
0161     // loop over ietas
0162     for (unsigned short absIeta = 1; absIeta <= 32; absIeta++) {
0163       hcal1 = caloTPGTranscoder->hcaletValue(absIeta, input);   // no eta-
0164       hcal2 = caloHcalScale->et(input, absIeta, 1);             // sign in transcoder
0165       hcal3 = caloTPGTranscoder->hcaletValue(-absIeta, input);  // no eta-
0166       hcal4 = caloHcalScale->et(input, absIeta, -1);            // sign in transcoder
0167 
0168       if ((!(hcal1 == hcal2)) || (!(hcal3 == hcal4))) {
0169         hcalIsConsistent = false;
0170         /*LogWarning("InconsistentData") 
0171          << "HCAL scales not consistent!"
0172          << "old HCAL is " << hcal1
0173          << " new HCAL is " << hcal2 ; */
0174       }
0175     }
0176   }
0177   if (!hcalIsConsistent) {
0178     // do something
0179     //cout << "WARNING: HCAL scales not consistent!" << endl;
0180     LogWarning("InconsistentData") << "HCAL scales not consistent!";
0181   } else {
0182     // do something else
0183     //cout << "HCAL scales okay" << endl;
0184   }
0185 }
0186 
0187 // ------------ method called once each job just before starting event loop  ------------
0188 void L1CaloInputScaleTester::beginJob() {}
0189 
0190 // ------------ method called once each job just after ending the event loop  ------------
0191 void L1CaloInputScaleTester::endJob() {}