CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/src/DQM/RCTMonitor/src/RCTMonitor.cc

Go to the documentation of this file.
00001 #include "DQM/RCTMonitor/interface/RCTMonitor.h"
00002 #include "DQM/RCTMonitor/interface/somedefinitions.h"
00003 #include "DQMServices/Core/interface/DQMStore.h"
00004 #include <iostream>
00005 
00006 RCTMonitor::RCTMonitor( const edm::ParameterSet& iConfig ):
00007   m_nevts(0),
00008   m_dbe(edm::Service<DQMStore>().operator->()),
00009   m_enableMonitorDaemon(iConfig.getUntrackedParameter<bool>("EnableMonitorDaemon")),
00010   m_rctSource(iConfig.getUntrackedParameter<edm::InputTag>("rctSource")),
00011   m_writeOutputFile(iConfig.getUntrackedParameter<bool>("WriteOutputFile")),
00012   m_outputFileName(iConfig.getUntrackedParameter<std::string>("OutputFileName"))
00013 {
00014 }
00015 
00016 
00017 RCTMonitor::~RCTMonitor()
00018 {
00019 }
00020 
00021 
00022 void RCTMonitor::beginJob()
00023 {
00024    BookRCT() ;
00025    
00026 }
00027 
00028 
00029 void RCTMonitor::endJob(void)
00030 {
00031   // Print out directory structure
00032   m_dbe->showDirStructure();
00033 
00034   // If requested write output to a root file
00035   if (m_writeOutputFile){
00036     m_dbe->save(m_outputFileName);
00037   }
00038 }
00039 
00040 
00041 void RCTMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup )
00042 {
00043 
00044   // Fill histograms
00045 
00046  
00047    FillRCT(iEvent,iSetup) ;
00048 
00049    
00050   // Increment number of events
00051   m_nevts++;
00052 
00053 }
00054 
00055 
00056 float DynamicScale(int EtaStamp)
00057 {
00058    //This function weights bin elements according to spatial extent of calorimeter tower.
00059    if(EtaStamp >= 6 && EtaStamp <= 15) {return ScaleINNER;}
00060    else if(EtaStamp==5 || EtaStamp==16) {return ScaleIN;}
00061    else if(EtaStamp == 4 || EtaStamp == 17) {return ScaleOUT;}
00062    else {return 0.000000;}
00063 }
00064 
00065 
00066    void RCTMonitor::FillRCT(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00067 {
00068 
00069   // Get the RCT digis
00070   edm::Handle<L1CaloEmCollection> em;
00071  //  edm::Handle<L1CaloRegionCollection> rgn;
00072 
00073  // iEvent.getByType(em);
00074  // iEvent.getByType(rgn);
00075   
00076    iEvent.getByLabel(m_rctSource,em);
00077 
00078 
00079 
00080   // Regions
00081 //  for (L1CaloRegionCollection::const_iterator ireg=rgn->begin(); ireg!=rgn->end(); ireg++) {
00082 
00083 //  if(ireg->et()>7){
00084 //    m_rctRegionsOccEtaPhi->Fill(ireg->gctPhi(),ireg->gctEta(),DynamicScale(ireg->gctEta()));
00085 //    m_rctRegionsEtEtaPhi->Fill(ireg->gctPhi(),ireg->gctEta(),ireg->et());
00086 //    m_rctRegionEt->Fill(ireg->et());
00087 //    m_rctTauVetoEtaPhi->Fill(ireg->gctPhi(),ireg->gctEta(),ireg->tauVeto());
00088 //   }
00089 //  }
00090 
00091 
00092   //Isolated and non-isolated EM with cut at >1 GeV
00093   for (L1CaloEmCollection::const_iterator iem=em->begin(); iem!=em->end(); iem++) {
00094    if(iem->rank()>1.){  //applies the 1 GeV cut
00095     if (iem->isolated()){  //looks for isolated EM candidates only
00096       m_rctIsoEmRank1->Fill(iem->rank());
00097 //std::cout << "Just to show what is there " << iem->rank() <<  std::endl ;
00098       m_rctIsoEmRankEtaPhi1->Fill(iem->regionId().iphi(),iem->regionId().ieta(),iem->rank());
00099       m_rctIsoEmOccEtaPhi1->Fill(iem->regionId().iphi(),iem->regionId().ieta(),DynamicScale(iem->regionId().ieta()));
00100       m_rctRelaxedEmRankEtaPhi1->Fill(iem->regionId().iphi(),iem->regionId().ieta(),iem->rank());
00101       m_rctRelaxedEmOccEtaPhi1->Fill(iem->regionId().iphi(),iem->regionId().ieta(),DynamicScale(iem->regionId().ieta()));
00102       m_rctRelaxedEmRank1->Fill(iem->rank());
00103     } else {  //instructions for Non-isolated EM candidates
00104       m_rctNonIsoEmRank1->Fill(iem->rank());
00105       m_rctNonIsoEmRankEtaPhi1->Fill(iem->regionId().iphi(),iem->regionId().ieta(),iem->rank());
00106       m_rctNonIsoEmOccEtaPhi1->Fill(iem->regionId().iphi(),iem->regionId().ieta(),DynamicScale(iem->regionId().ieta()));
00107       m_rctRelaxedEmRankEtaPhi1->Fill(iem->regionId().iphi(),iem->regionId().ieta(),iem->rank());
00108       m_rctRelaxedEmOccEtaPhi1->Fill(iem->regionId().iphi(),iem->regionId().ieta(),DynamicScale(iem->regionId().ieta()));
00109       m_rctRelaxedEmRank1->Fill(iem->rank());
00110     }
00111    }
00112    if(iem->rank()>10.){  //applies the 10 GeV cut
00113     if (iem->isolated()){  //looks for isolated EM candidates only
00114       m_rctIsoEmOccEtaPhi10->Fill(iem->regionId().iphi(),iem->regionId().ieta(),DynamicScale(iem->regionId().ieta()));
00115       m_rctRelaxedEmOccEtaPhi10->Fill(iem->regionId().iphi(),iem->regionId().ieta(),DynamicScale(iem->regionId().ieta()));
00116     } else {  //instructions for Non-isolated EM candidates
00117       m_rctNonIsoEmOccEtaPhi10->Fill(iem->regionId().iphi(),iem->regionId().ieta(),DynamicScale(iem->regionId().ieta()));
00118       m_rctRelaxedEmOccEtaPhi10->Fill(iem->regionId().iphi(),iem->regionId().ieta(),DynamicScale(iem->regionId().ieta()));;
00119     }
00120    }
00121   }
00122 }
00123 
00124 
00125 
00126 void RCTMonitor::BookRCT()
00127 {
00128 //std::cout << "I am in the RCT booking"  << std::endl ;
00129 
00130   // Book RCT histograms
00131   m_dbe->setCurrentFolder("RCT");
00132 
00133   m_rctIsoEmRankEtaPhi1     = m_dbe->book2D("RctIsoEmRankEtaPhi",      "ISO EM RANK"         , PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00134   m_rctIsoEmOccEtaPhi1      = m_dbe->book2D("RctIsoEmOccEtaPhi",       "ISO EM OCCUPANCY"    , PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00135   m_rctIsoEmRank1           = m_dbe->book1D("RctIsoEmRank",            "ISO EM RANK"         , R6BINS, R6MIN, R6MAX);
00136   m_rctIsoEmRankEtaPhi10    = m_dbe->book2D("RctIsoEmRankEtaPhi10",    "ISO EM RANK"         , PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00137   m_rctIsoEmOccEtaPhi10     = m_dbe->book2D("RctIsoEmOccEtaPhi10",     "ISO EM OCCUPANCY"    , PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00138   m_rctIsoEmRank10          = m_dbe->book1D("RctIsoEmRank10",          "ISO EM RANK"         , R6BINS, R6MIN, R6MAX);
00139 
00140   m_rctNonIsoEmRankEtaPhi1  = m_dbe->book2D("RctNonIsoEmRankEtaPhi",   "NON-ISO EM RANK"     , PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00141   m_rctNonIsoEmOccEtaPhi1   = m_dbe->book2D("RctNonIsoEmOccEtaPhi",    "NON-ISO EM OCCUPANCY", PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00142   m_rctNonIsoEmRank1        = m_dbe->book1D("RctNonIsoEmRank",         "NON-ISO EM RANK"     , R6BINS, R6MIN, R6MAX);
00143   m_rctNonIsoEmRankEtaPhi10 = m_dbe->book2D("RctNonIsoEmRankEtaPhi10", "NON-ISO EM RANK"     , PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00144   m_rctNonIsoEmOccEtaPhi10  = m_dbe->book2D("RctNonIsoEmOccEtaPhi10",  "NON-ISO EM OCCUPANCY", PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00145   m_rctNonIsoEmRank10       = m_dbe->book1D("RctNonIsoEmRank10",       "NON-ISO EM RANK"     , R6BINS, R6MIN, R6MAX);
00146 
00147   m_rctRelaxedEmRankEtaPhi1 = m_dbe->book2D("RctRelaxedEmRankEtaPhi",  "RELAXED EM RANK"     , PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00148   m_rctRelaxedEmOccEtaPhi1  = m_dbe->book2D("RctRelaxedEmOccEtaPhi",   "RELAXED EM OCCUPANCY", PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00149   m_rctRelaxedEmRank1       = m_dbe->book1D("RctRelaxedEmRank",        "RELAXED EM RANK"     , R6BINS, R6MIN, R6MAX);
00150   m_rctRelaxedEmRankEtaPhi10= m_dbe->book2D("RctRelaxedEmRankEtaPhi",  "RELAXED EM RANK"     , PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00151   m_rctRelaxedEmOccEtaPhi10 = m_dbe->book2D("RctRelaxedEmOccEtaPhi10",   "RELAXED EM OCCUPANCY", PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00152   m_rctRelaxedEmRank10      = m_dbe->book1D("RctRelaxedEmRank",        "RELAXED EM RANK"     , R6BINS, R6MIN, R6MAX);
00153 
00154   m_rctRegionsEtEtaPhi      = m_dbe->book2D("RctRegionsEtEtaPhi",      "REGION E_{T}"        , PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00155   m_rctRegionsOccEtaPhi     = m_dbe->book2D("RctRegionsOccEtaPhi",     "REGION OCCUPANCY"    , PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00156   m_rctTauVetoEtaPhi        = m_dbe->book2D("RctTauVetoEtaPhi",        "TAU VETO OCCUPANCY"  , PHIBINS, PHIMIN, PHIMAX, ETABINS, ETAMIN, ETAMAX);
00157   m_rctRegionEt             = m_dbe->book1D("RctRegionEt",             "REGION E_{T}"        , R10BINS, R10MIN, R10MAX);
00158 
00159 
00160 
00161 }
00162 
00163 
00164 // define this as a plug-in
00165 DEFINE_FWK_MODULE(RCTMonitor);
00166 
00167   
00168