CMS 3D CMS Logo

HLTMuonL1DQMSource.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    HLTMuonL1DQMSource
00004 // Class:      HLTMuonL1DQMSource
00005 // 
00013 //
00014 // Original Author:  Muriel VANDER DONCKT *:0
00015 //         Created:  Wed Dec 12 09:55:42 CET 2007
00016 // $Id: HLTMuonL1DQMSource.cc,v 1.1 2008/06/25 10:46:57 muriel Exp $
00017 //
00018 //
00019 
00020 
00021 
00022 #include "DQM/HLTEvF/interface/HLTMuonL1DQMSource.h"
00023 #include "DQMServices/Core/interface/DQMStore.h"
00024 #include "DataFormats/Common/interface/Handle.h"
00025 #include "FWCore/Framework/interface/Event.h"
00026 #include "FWCore/Framework/interface/EventSetup.h"
00027 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00028 #include "FWCore/ServiceRegistry/interface/Service.h"
00029 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
00030 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
00031 #include "DataFormats/TrackReco/interface/Track.h"
00032 #include "DataFormats/L1Trigger/interface/L1MuonParticle.h"
00033 #include "DataFormats/L1Trigger/interface/L1MuonParticleFwd.h"
00034 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTCand.h"
00035 
00036 #include "TMath.h" 
00037 
00038 
00039 using namespace std;
00040 using namespace edm;
00041 using namespace l1extra;
00042 //
00043 // constructors and destructor
00044 //
00045 HLTMuonL1DQMSource::HLTMuonL1DQMSource( const edm::ParameterSet& parameters_ ) :counterEvt_(0)
00046 
00047 {
00048   verbose_ = parameters_.getUntrackedParameter < bool > ("verbose", false);
00049   monitorName_ = parameters_.getUntrackedParameter<string>("monitorName","HLT/HLTMuon");
00050   prescaleEvt_ = parameters_.getUntrackedParameter<int>("prescaleEvt", -1);
00051   level_ = parameters_.getUntrackedParameter<int>("Level",2);
00052   l1muCollectionTag_ = parameters_.getUntrackedParameter<InputTag>("L1MuonTag",edm::InputTag("hltL1extraParticles"));
00053 
00054    dbe_ = 0 ;
00055    if (parameters_.getUntrackedParameter < bool > ("DQMStore", false)) {
00056      dbe_ = Service < DQMStore > ().operator->();
00057      dbe_->setVerbose(0);
00058    }
00059  
00060    outputFile_ =
00061        parameters_.getUntrackedParameter < std::string > ("outputFile", "");
00062    if (outputFile_.size() != 0 && verbose_) {
00063      std::cout << "Muon HLT Monitoring histograms will be saved to " 
00064                << outputFile_ << std::endl;
00065    }
00066    else {
00067      outputFile_ = "HLTMuonDQM.root";
00068    }
00069  
00070    bool disable =
00071      parameters_.getUntrackedParameter < bool > ("disableROOToutput", false);
00072    if (disable) {
00073      outputFile_ = "";
00074    }
00075  
00076    if (dbe_ != NULL) {
00077      dbe_->setCurrentFolder("HLT/HLTMonMuon");
00078    }
00079 
00080 
00081 }
00082 
00083 
00084 HLTMuonL1DQMSource::~HLTMuonL1DQMSource()
00085 {
00086    
00087   // do anything here that needs to be done at desctruction time
00088   // (e.g. close files, deallocate resources etc.)
00089   
00090 }
00091 
00092 //--------------------------------------------------------
00093 void HLTMuonL1DQMSource::beginJob(const EventSetup& context){
00094 
00095  
00096    if (dbe_) {
00097      dbe_->setCurrentFolder(monitorName_);
00098      if (monitorName_ != "" ) monitorName_ = monitorName_+"/" ;
00099      if (verbose_)cout << "===>DQM event prescale = " << prescaleEvt_ << " events "<< endl;
00100      
00101      
00103      const int NBINS = 100; 
00104 
00105      // create and cd into new folder
00106      char name[512], title[512];
00107      sprintf(name,"Level%i",level_);
00108      dbe_->setCurrentFolder(monitorName_+name);
00109      hl1quality = dbe_->book1D("h1L1Quality","GMT quality Flag", 8, 0., 8.);
00110      sprintf(name,"HLTMuonL%i_NMu",level_);
00111      sprintf(title,"L%i number of muons",level_);
00112      hNMu = dbe_->book1D(name,title, 5, 0., 5.);
00113      hNMu->setAxisTitle("Number of muons", 1);
00114      sprintf(name,"HLTMuonL%i_pt",level_);
00115      sprintf(title,"L%i Pt",level_);
00116      hpt = dbe_->book1D(name,title, NBINS, 0., 100);
00117      hpt->setAxisTitle("Pt", 1);
00118      sprintf(name,"HLTMuonL%i_eta",level_);
00119      sprintf(title,"L%i Muon #eta",level_);
00120      heta = dbe_->book1D(name,title, NBINS, -2.5, 2.5);
00121      heta->setAxisTitle("#eta", 1);
00122      sprintf(name,"HLTMuonL%i_phi",level_);
00123      sprintf(title,"L%i Muon #phi",level_);
00124      hphi = dbe_->book1D(name,title, NBINS, -3.15, 3.15);
00125      hphi->setAxisTitle("#phi", 1);
00126      sprintf(name,"HLTMuonL%i_etaphi",level_);
00127      sprintf(title,"L%i Muon #eta vs #phi",level_);
00128      hetaphi = dbe_->book2D(name,title, NBINS, -3.15, 3.15,NBINS,-2.5, 2.5);
00129      hetaphi->setAxisTitle("#phi", 1);
00130      hetaphi->setAxisTitle("#eta", 2); 
00131      sprintf(name,"HLTMuonL%i_ptphi",level_);
00132      sprintf(title,"L%i Muon pt vs #phi",level_);         
00133      hptphi = dbe_->book2D(name,title, NBINS, 0., 100.,NBINS,-3.15, 3.15);
00134      hptphi->setAxisTitle("pt", 1);
00135      hptphi->setAxisTitle("#phi", 2);
00136      sprintf(name,"HLTMuonL%i_pteta",level_);
00137      sprintf(title,"L%i Muon pt vs #eta",level_);         
00138      hpteta = dbe_->book2D(name,title, NBINS, 0., 100.,NBINS,-2.5, 2.5);
00139      hpteta->setAxisTitle("pt", 1);
00140      hpteta->setAxisTitle("#eta", 2);
00141      sprintf(name,"HLTMuonL%i_charge",level_);
00142      sprintf(title,"L%i Muon Charge",level_);         
00143      hcharge  = dbe_->book1D(name,title, 3, -1.5, 1.5);
00144      hcharge->setAxisTitle("Charge", 1);
00145      if(verbose_)dbe_->showDirStructure();
00146        
00147      // Muon det id is 2 pushed in bits 28:31
00148      const unsigned int detector_id = 2<<28;
00149      dbe_->tagContents(monitorName_, detector_id);
00150    } 
00151 }
00152 
00153 //--------------------------------------------------------
00154 void HLTMuonL1DQMSource::beginRun(const edm::Run& r, const EventSetup& context) {
00155   // reset all me's
00156   vector<MonitorElement*> AllME=dbe_->getAllContents(monitorName_);
00157   vector<MonitorElement*>::iterator me=AllME.begin();
00158   for ( ; me != AllME.end() ; ++me ){
00159     (*me)->Reset();
00160   }
00161 }
00162 
00163 //--------------------------------------------------------
00164 void HLTMuonL1DQMSource::beginLuminosityBlock(const LuminosityBlock& lumiSeg, 
00165                                       const EventSetup& context) {
00166   
00167 }
00168 
00169 // ----------------------------------------------------------
00170 void HLTMuonL1DQMSource::analyze(const Event& iEvent, 
00171                          const EventSetup& iSetup )
00172 {  
00173   if ( !dbe_) return;
00174   counterEvt_++;
00175   if (prescaleEvt_ > 0 && counterEvt_%prescaleEvt_!=0) return;
00176   if (verbose_)cout << " processing conterEvt_: " << counterEvt_ <<endl;
00177 
00178   edm::Handle<L1MuonParticleCollection> muColl;
00179   iEvent.getByLabel(l1muCollectionTag_, muColl);
00180   if (!muColl.failedToGet()){
00181     L1MuonParticleCollection::const_iterator l1ref;
00182     L1MuonParticleRef::key_type l1ParticleIndex = 0;
00183     hNMu->Fill(muColl->size());
00184     
00185     for(l1ref = muColl->begin(); l1ref != muColl->end(); ++l1ref,++l1ParticleIndex) {    
00186       hcharge->Fill(l1ref->charge()); 
00187       hpt->Fill(l1ref->pt());
00188       hphi->Fill(l1ref->phi());
00189       heta->Fill(l1ref->eta());
00190       hetaphi->Fill(l1ref->phi(),l1ref->eta());
00191       hptphi->Fill(l1ref->pt(),l1ref->phi());
00192       hpteta->Fill(l1ref->pt(),l1ref->eta());
00193       hl1quality->Fill(l1ref->gmtMuonCand().quality());
00194     }
00195   }
00196 }
00197 
00198 
00199 
00200 //--------------------------------------------------------
00201 void HLTMuonL1DQMSource::endLuminosityBlock(const LuminosityBlock& lumiSeg, 
00202                                     const EventSetup& context) {
00203 }
00204 //--------------------------------------------------------
00205 void HLTMuonL1DQMSource::endRun(const Run& r, const EventSetup& context){
00206 }
00207 //--------------------------------------------------------
00208 void HLTMuonL1DQMSource::endJob(){
00209    LogInfo("HLTMonMuon") << "analyzed " << counterEvt_ << " events";
00210  
00211    if (outputFile_.size() != 0 && dbe_)
00212     dbe_->save(outputFile_);
00213  
00214    return;
00215 }

Generated on Tue Jun 9 17:33:08 2009 for CMSSW by  doxygen 1.5.4