CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/L1TriggerConfig/L1ScalesProducers/src/L1CaloInputScaleTester.cc

Go to the documentation of this file.
00001 #include "L1TriggerConfig/L1ScalesProducers/interface/L1CaloInputScaleTester.h"
00002 
00003 // system include files
00004 #include <memory>
00005 #include <iostream>
00006 using std::cout;
00007 using std::endl;
00008 #include <iomanip>
00009 using std::setprecision;
00010 
00011 // user include files
00012 #include "FWCore/Framework/interface/Frameworkfwd.h"
00013 
00014 #include "FWCore/Framework/interface/MakerMacros.h"
00015 #include "FWCore/Framework/interface/EventSetup.h"
00016 #include "FWCore/Framework/interface/ESHandle.h"
00017 
00018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00019 
00020 #include "CalibFormats/CaloTPG/interface/CaloTPGTranscoder.h"
00021 #include "CalibFormats/CaloTPG/interface/CaloTPGRecord.h"
00022 #include "CalibCalorimetry/EcalTPGTools/interface/EcalTPGScale.h"
00023 #include "CondFormats/L1TObjects/interface/L1CaloEcalScale.h"
00024 #include "CondFormats/DataRecord/interface/L1CaloEcalScaleRcd.h"
00025 #include "CondFormats/L1TObjects/interface/L1CaloHcalScale.h"
00026 #include "CondFormats/DataRecord/interface/L1CaloHcalScaleRcd.h"
00027 #include "DataFormats/EcalDetId/interface/EcalTrigTowerDetId.h"
00028 
00029 
00030 //
00031 // constants, enums and typedefs
00032 //
00033 
00034 //
00035 // static data member definitions
00036 //
00037 
00038 //
00039 // constructors and destructor
00040 //
00041 L1CaloInputScaleTester::L1CaloInputScaleTester(const edm::ParameterSet& iConfig)
00042 
00043 {
00044    //now do what ever initialization is needed
00045 
00046 }
00047 
00048 
00049 L1CaloInputScaleTester::~L1CaloInputScaleTester()
00050 {
00051  
00052    // do anything here that needs to be done at destruction time
00053    // (e.g. close files, deallocate resources etc.)
00054 
00055 }
00056 
00057 
00058 //
00059 // member functions
00060 //
00061 
00062 // ------------ method called to for each event  ------------
00063 void
00064 L1CaloInputScaleTester::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00065 {
00066    using namespace edm;
00067 
00068 #ifdef THIS_IS_AN_EVENT_EXAMPLE
00069    Handle<ExampleData> pIn;
00070    iEvent.getByLabel("example",pIn);
00071 #endif
00072    
00073 #ifdef THIS_IS_AN_EVENTSETUP_EXAMPLE
00074    ESHandle<SetupData> pSetup;
00075    iSetup.get<SetupRecord>().get(pSetup);
00076 #endif
00077 
00078    ESHandle<L1CaloEcalScale> caloEcalScale;
00079    ESHandle<L1CaloHcalScale> caloHcalScale;
00080    ESHandle<CaloTPGTranscoder> caloTPGTranscoder;
00081    iSetup.get<L1CaloEcalScaleRcd>().get(caloEcalScale);
00082    iSetup.get<L1CaloHcalScaleRcd>().get(caloHcalScale);
00083    iSetup.get<CaloTPGRecord>().get(caloTPGTranscoder);
00084    
00085    EcalTPGScale* ecalTPGScale = new EcalTPGScale();
00086    ecalTPGScale->setEventSetup(iSetup);
00087 
00088    bool ecalIsConsistent = true;
00089    bool hcalIsConsistent = true;
00090 
00091    double ecal1; 
00092    double ecal2;
00093    double hcal1;
00094    double hcal2;
00095    double hcal3;
00096    double hcal4;
00097 
00098 
00099    // compare the ecal scales
00100 
00101    // 8 bits of input energy
00102    for (unsigned short input = 0; input <= 0xFF; input++)
00103      {
00104        // loop over ietas, barrel first
00105        for (unsigned short absIeta = 1; absIeta <= 17; absIeta++)
00106          {
00107 
00108            // positive eta
00109            ecal1 = ecalTPGScale->getTPGInGeV( (unsigned int) input, 
00110                                               EcalTrigTowerDetId(1, EcalBarrel,
00111                                                                  absIeta, 1));
00112            ecal2 = caloEcalScale->et(input, absIeta, 1);
00113 
00114            if ( !(ecal1 == ecal2) )
00115              {
00116                ecalIsConsistent = false;
00117                /*LogWarning("InconsistentData") 
00118                  << "ECAL scales not consistent! (pos eta barrel)"
00119                  << "old ECAL is " << setprecision (8) << ecal1
00120                  << " new ECAL is " << setprecision (8) << ecal2 
00121                  << " difference is " << (ecal1 - ecal2) ; */
00122              }
00123 
00124            // negative eta
00125            ecal1 = ecalTPGScale->
00126              getTPGInGeV( (unsigned int) input, EcalTrigTowerDetId(-1, EcalBarrel,
00127                                                            absIeta, 2));
00128            ecal2 = caloEcalScale->et(input, absIeta, -1);
00129 
00130            if ( !(ecal1 == ecal2) )
00131              {
00132                ecalIsConsistent = false;
00133                /*LogWarning("InconsistentData") 
00134                  << "ECAL scales not consistent! (neg eta barrel)"
00135                  << "old ECAL is " << setprecision (8) << ecal1
00136                  << " new ECAL is " << setprecision (8) << ecal2               
00137                  << " difference is " << (ecal1 - ecal2) ; */       
00138              }
00139          }
00140        // now loop over endcap ietas
00141        for (unsigned short absIeta = 18; absIeta <= 28; absIeta++)
00142          {
00143 
00144            // positive eta
00145            ecal1 = ecalTPGScale->
00146              getTPGInGeV( (unsigned int) input, EcalTrigTowerDetId(1, EcalEndcap,
00147                                                            absIeta, 1));
00148            ecal2 = caloEcalScale->et(input, absIeta, 1);
00149 
00150            if ( !(ecal1 == ecal2) )
00151              {
00152                ecalIsConsistent = false;
00153                /*LogWarning("InconsistentData") 
00154                  << "ECAL scales not consistent! (pos eta endcap)"
00155                  << "old ECAL is " << setprecision (8) << ecal1
00156                  << " new ECAL is " << setprecision (8) << ecal2 
00157                  << " difference is " << (ecal1 - ecal2) ; */
00158              }
00159            
00160            // negative eta
00161            ecal1 = ecalTPGScale->
00162              getTPGInGeV( (unsigned int) input, EcalTrigTowerDetId(-1, EcalEndcap,
00163                                                            absIeta, 2));
00164            ecal2 = caloEcalScale->et(input, absIeta, -1);
00165 
00166            if ( !(ecal1 == ecal2) )
00167              {
00168                ecalIsConsistent = false;
00169                /*LogWarning("InconsistentData") 
00170                  << "ECAL scales not consistent! (neg eta endcap)"
00171                  << "old ECAL is " << setprecision (8) << ecal1
00172                  << " new ECAL is " << setprecision (8) << ecal2 
00173                  << " difference is " << (ecal1 - ecal2) ; */
00174              }
00175          }
00176      }
00177 
00178    if (!ecalIsConsistent)
00179      {
00180        // do something
00181        //cout << "WARNING: ECAL scales not consistent!" << endl;
00182        LogWarning("InconsistentData") << "ECAL scales not consistent!";
00183      }
00184    else
00185      {
00186        // do something else
00187        //cout << "ECAL scales okay" << endl;
00188      }
00189 
00190    // compare the hcal scales
00191 
00192    for (unsigned short input = 0; input <= 0xFF; input++)
00193      {
00194        // loop over ietas
00195        for (unsigned short absIeta = 1; absIeta <= 32; absIeta++)
00196          {
00197            hcal1 = caloTPGTranscoder->hcaletValue(absIeta, input); // no eta-
00198            hcal2 = caloHcalScale->et(input, absIeta, 1); // sign in transcoder
00199            hcal3 = caloTPGTranscoder->hcaletValue(-absIeta, input); // no eta-
00200            hcal4 = caloHcalScale->et(input, absIeta,-1); // sign in transcoder
00201 
00202 
00203 
00204            if ( (!(hcal1 == hcal2))||(!(hcal3==hcal4)))
00205              {
00206                hcalIsConsistent = false;
00207                /*LogWarning("InconsistentData") 
00208                  << "HCAL scales not consistent!"
00209                  << "old HCAL is " << hcal1
00210                  << " new HCAL is " << hcal2 ; */
00211              }
00212          }
00213      }
00214    if (!hcalIsConsistent)
00215      {
00216        // do something
00217        //cout << "WARNING: HCAL scales not consistent!" << endl;
00218        LogWarning("InconsistentData") << "HCAL scales not consistent!";
00219      }
00220    else
00221      {
00222        // do something else
00223        //cout << "HCAL scales okay" << endl;
00224      }
00225 }
00226 
00227 
00228 // ------------ method called once each job just before starting event loop  ------------
00229 void 
00230 L1CaloInputScaleTester::beginJob()
00231 {
00232 }
00233 
00234 // ------------ method called once each job just after ending the event loop  ------------
00235 void 
00236 L1CaloInputScaleTester::endJob() {
00237 }
00238