CMS 3D CMS Logo

Public Member Functions | Private Member Functions

L1CaloInputScaleTester Class Reference

#include <L1TriggerConfig/L1ScalesProducers/src/L1CaloInputScaleTester.cc>

Inheritance diagram for L1CaloInputScaleTester:
edm::EDAnalyzer edm::EDConsumerBase

List of all members.

Public Member Functions

 L1CaloInputScaleTester (const edm::ParameterSet &)
 ~L1CaloInputScaleTester ()

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob ()
virtual void endJob ()

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 28 of file L1CaloInputScaleTester.h.


Constructor & Destructor Documentation

L1CaloInputScaleTester::L1CaloInputScaleTester ( const edm::ParameterSet iConfig) [explicit]

Definition at line 41 of file L1CaloInputScaleTester.cc.

{
   //now do what ever initialization is needed

}
L1CaloInputScaleTester::~L1CaloInputScaleTester ( )

Definition at line 49 of file L1CaloInputScaleTester.cc.

{
 
   // do anything here that needs to be done at destruction time
   // (e.g. close files, deallocate resources etc.)

}

Member Function Documentation

void L1CaloInputScaleTester::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [private, virtual]

Implements edm::EDAnalyzer.

Definition at line 64 of file L1CaloInputScaleTester.cc.

References EcalBarrel, EcalEndcap, edm::EventSetup::get(), edm::Event::getByLabel(), and LaserDQM_cfg::input.

{
   using namespace edm;

#ifdef THIS_IS_AN_EVENT_EXAMPLE
   Handle<ExampleData> pIn;
   iEvent.getByLabel("example",pIn);
#endif
   
#ifdef THIS_IS_AN_EVENTSETUP_EXAMPLE
   ESHandle<SetupData> pSetup;
   iSetup.get<SetupRecord>().get(pSetup);
#endif

   ESHandle<L1CaloEcalScale> caloEcalScale;
   ESHandle<L1CaloHcalScale> caloHcalScale;
   ESHandle<CaloTPGTranscoder> caloTPGTranscoder;
   iSetup.get<L1CaloEcalScaleRcd>().get(caloEcalScale);
   iSetup.get<L1CaloHcalScaleRcd>().get(caloHcalScale);
   iSetup.get<CaloTPGRecord>().get(caloTPGTranscoder);
   
   EcalTPGScale* ecalTPGScale = new EcalTPGScale();
   ecalTPGScale->setEventSetup(iSetup);

   bool ecalIsConsistent = true;
   bool hcalIsConsistent = true;

   double ecal1; 
   double ecal2;
   double hcal1;
   double hcal2;
   double hcal3;
   double hcal4;


   // compare the ecal scales

   // 8 bits of input energy
   for (unsigned short input = 0; input <= 0xFF; input++)
     {
       // loop over ietas, barrel first
       for (unsigned short absIeta = 1; absIeta <= 17; absIeta++)
         {

           // positive eta
           ecal1 = ecalTPGScale->getTPGInGeV( (unsigned int) input, 
                                              EcalTrigTowerDetId(1, EcalBarrel,
                                                                 absIeta, 1));
           ecal2 = caloEcalScale->et(input, absIeta, 1);

           if ( !(ecal1 == ecal2) )
             {
               ecalIsConsistent = false;
               /*LogWarning("InconsistentData") 
                 << "ECAL scales not consistent! (pos eta barrel)"
                 << "old ECAL is " << setprecision (8) << ecal1
                 << " new ECAL is " << setprecision (8) << ecal2 
                 << " difference is " << (ecal1 - ecal2) ; */
             }

           // negative eta
           ecal1 = ecalTPGScale->
             getTPGInGeV( (unsigned int) input, EcalTrigTowerDetId(-1, EcalBarrel,
                                                           absIeta, 2));
           ecal2 = caloEcalScale->et(input, absIeta, -1);

           if ( !(ecal1 == ecal2) )
             {
               ecalIsConsistent = false;
               /*LogWarning("InconsistentData") 
                 << "ECAL scales not consistent! (neg eta barrel)"
                 << "old ECAL is " << setprecision (8) << ecal1
                 << " new ECAL is " << setprecision (8) << ecal2               
                 << " difference is " << (ecal1 - ecal2) ; */       
             }
         }
       // now loop over endcap ietas
       for (unsigned short absIeta = 18; absIeta <= 28; absIeta++)
         {

           // positive eta
           ecal1 = ecalTPGScale->
             getTPGInGeV( (unsigned int) input, EcalTrigTowerDetId(1, EcalEndcap,
                                                           absIeta, 1));
           ecal2 = caloEcalScale->et(input, absIeta, 1);

           if ( !(ecal1 == ecal2) )
             {
               ecalIsConsistent = false;
               /*LogWarning("InconsistentData") 
                 << "ECAL scales not consistent! (pos eta endcap)"
                 << "old ECAL is " << setprecision (8) << ecal1
                 << " new ECAL is " << setprecision (8) << ecal2 
                 << " difference is " << (ecal1 - ecal2) ; */
             }
           
           // negative eta
           ecal1 = ecalTPGScale->
             getTPGInGeV( (unsigned int) input, EcalTrigTowerDetId(-1, EcalEndcap,
                                                           absIeta, 2));
           ecal2 = caloEcalScale->et(input, absIeta, -1);

           if ( !(ecal1 == ecal2) )
             {
               ecalIsConsistent = false;
               /*LogWarning("InconsistentData") 
                 << "ECAL scales not consistent! (neg eta endcap)"
                 << "old ECAL is " << setprecision (8) << ecal1
                 << " new ECAL is " << setprecision (8) << ecal2 
                 << " difference is " << (ecal1 - ecal2) ; */
             }
         }
     }

   if (!ecalIsConsistent)
     {
       // do something
       //cout << "WARNING: ECAL scales not consistent!" << endl;
       LogWarning("InconsistentData") << "ECAL scales not consistent!";
     }
   else
     {
       // do something else
       //cout << "ECAL scales okay" << endl;
     }

   // compare the hcal scales

   for (unsigned short input = 0; input <= 0xFF; input++)
     {
       // loop over ietas
       for (unsigned short absIeta = 1; absIeta <= 32; absIeta++)
         {
           hcal1 = caloTPGTranscoder->hcaletValue(absIeta, input); // no eta-
           hcal2 = caloHcalScale->et(input, absIeta, 1); // sign in transcoder
           hcal3 = caloTPGTranscoder->hcaletValue(-absIeta, input); // no eta-
           hcal4 = caloHcalScale->et(input, absIeta,-1); // sign in transcoder



           if ( (!(hcal1 == hcal2))||(!(hcal3==hcal4)))
             {
               hcalIsConsistent = false;
               /*LogWarning("InconsistentData") 
                 << "HCAL scales not consistent!"
                 << "old HCAL is " << hcal1
                 << " new HCAL is " << hcal2 ; */
             }
         }
     }
   if (!hcalIsConsistent)
     {
       // do something
       //cout << "WARNING: HCAL scales not consistent!" << endl;
       LogWarning("InconsistentData") << "HCAL scales not consistent!";
     }
   else
     {
       // do something else
       //cout << "HCAL scales okay" << endl;
     }
}
void L1CaloInputScaleTester::beginJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 230 of file L1CaloInputScaleTester.cc.

{
}
void L1CaloInputScaleTester::endJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 236 of file L1CaloInputScaleTester.cc.

                               {
}