CMS 3D CMS Logo

L1CaloInputScalesGenerator.cc

Go to the documentation of this file.
00001 #include "L1TriggerConfig/L1ScalesProducers/interface/L1CaloInputScalesGenerator.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 #include <fstream>
00011 using std::ofstream;
00012 
00013 
00014 // user include files
00015 #include "FWCore/Framework/interface/Frameworkfwd.h"
00016 
00017 #include "FWCore/Framework/interface/MakerMacros.h"
00018 #include "FWCore/Framework/interface/EventSetup.h"
00019 #include "FWCore/Framework/interface/ESHandle.h"
00020 
00021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00022 
00023 #include "CalibFormats/CaloTPG/interface/CaloTPGTranscoder.h"
00024 #include "CalibFormats/CaloTPG/interface/CaloTPGRecord.h"
00025 #include "CalibCalorimetry/EcalTPGTools/interface/EcalTPGScale.h"
00026 #include "DataFormats/EcalDetId/interface/EcalTrigTowerDetId.h"
00027 
00028 
00029 //
00030 // constants, enums and typedefs
00031 //
00032 
00033 //
00034 // static data member definitions
00035 //
00036 
00037 //
00038 // constructors and destructor
00039 //
00040 L1CaloInputScalesGenerator::L1CaloInputScalesGenerator(const edm::ParameterSet& iConfig)
00041 
00042 {
00043    //now do what ever initialization is needed
00044 
00045 }
00046 
00047 
00048 L1CaloInputScalesGenerator::~L1CaloInputScalesGenerator()
00049 {
00050  
00051    // do anything here that needs to be done at desctruction time
00052    // (e.g. close files, deallocate resources etc.)
00053 
00054 }
00055 
00056 
00057 //
00058 // member functions
00059 //
00060 
00061 // ------------ method called to for each event  ------------
00062 void
00063 L1CaloInputScalesGenerator::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00064 {
00065    using namespace edm;
00066 
00067 #ifdef THIS_IS_AN_EVENT_EXAMPLE
00068    Handle<ExampleData> pIn;
00069    iEvent.getByLabel("example",pIn);
00070 #endif
00071    
00072 #ifdef THIS_IS_AN_EVENTSETUP_EXAMPLE
00073    ESHandle<SetupData> pSetup;
00074    iSetup.get<SetupRecord>().get(pSetup);
00075 #endif
00076 
00077    ESHandle<CaloTPGTranscoder> caloTPGTranscoder;
00078    iSetup.get<CaloTPGRecord>().get(caloTPGTranscoder);
00079    
00080    EcalTPGScale* ecalTPGScale = new EcalTPGScale();
00081    ecalTPGScale->setEventSetup(iSetup);
00082 
00083    double output; 
00084    ofstream scalesFile("L1CaloInputScales.cfi");
00085 
00086 
00087    // Write the ecal scales, positive eta
00088 
00089    scalesFile << "es_module L1CaloInputScalesProducer = "
00090               << "L1CaloInputScalesProducer {\n\tvdouble "
00091               << "L1EcalEtThresholdsPositiveEta = {" << endl; 
00092 
00093    // loop over ietas, barrel
00094    for (unsigned short absIeta = 1; absIeta <= 28; absIeta++)
00095      {
00096        EcalSubdetector subdet = ( absIeta <= 17 ) ? EcalBarrel : EcalEndcap ;
00097        // 8 bits of input energy
00098        for (unsigned short input = 0; input <= 0xFF; input++)
00099          {
00100            output = ecalTPGScale->getTPGInGeV( (uint) input, 
00101                                               EcalTrigTowerDetId(1, subdet,
00102                                                                  absIeta, 1));
00103            scalesFile << setprecision (8) << output;
00104            if (absIeta == 28 && input == 0xFF)
00105              {
00106                scalesFile << "}";
00107              }
00108            else
00109              {
00110                scalesFile << ", ";
00111              }
00112          }
00113        scalesFile << endl;
00114      }
00115 
00116    // Write the ecal scales, negative eta
00117 
00118    scalesFile << endl << "\tvdouble L1EcalEtThresholdsNegativeEta = {" << endl;
00119 
00120    // loop over ietas, barrel first
00121    for (unsigned short absIeta = 1; absIeta <= 28; absIeta++)
00122      {
00123        EcalSubdetector subdet = ( absIeta <= 17 ) ? EcalBarrel : EcalEndcap ;
00124        // 8 bits of input energy
00125        for (unsigned short input = 0; input <= 0xFF; input++)
00126          {
00127            // negative eta
00128            output = ecalTPGScale->
00129              getTPGInGeV( (uint) input, EcalTrigTowerDetId(-1, subdet,
00130                                                            absIeta, 2));
00131            scalesFile << setprecision (8) << output;
00132            if (absIeta == 28 && input == 0xFF)
00133              {
00134                scalesFile << "}";
00135              }
00136            else
00137              {
00138                scalesFile << ", ";
00139              }
00140          }
00141        scalesFile << endl;
00142      }
00143 
00144    // Write the hcal scales
00145 
00146      scalesFile << endl << "\tvdouble L1HcalEtThresholds = {" << endl;
00147 
00148    // loop over ietas
00149    for (unsigned short absIeta = 1; absIeta <= 32; absIeta++)
00150      {
00151        for (unsigned short input = 0; input <= 0xFF; input++)
00152          {
00153            output = caloTPGTranscoder->hcaletValue(absIeta, input); 
00154            scalesFile << setprecision (8) << output ;
00155            if (absIeta == 32 && input == 0xFF)
00156              {
00157                scalesFile << "}";
00158              }
00159            else
00160              {
00161                scalesFile << ", ";
00162              }
00163          }
00164        scalesFile << endl;
00165      }
00166    scalesFile << "}" << endl;
00167 
00168    scalesFile.close();
00169 }
00170 
00171 
00172 // ------------ method called once each job just before starting event loop  ------------
00173 void 
00174 L1CaloInputScalesGenerator::beginJob(const edm::EventSetup&)
00175 {
00176 }
00177 
00178 // ------------ method called once each job just after ending the event loop  ------------
00179 void 
00180 L1CaloInputScalesGenerator::endJob() {
00181 }
00182 

Generated on Tue Jun 9 17:40:28 2009 for CMSSW by  doxygen 1.5.4