CMS 3D CMS Logo

HLTMuonDQMSource.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    HLTMuonDQMSource
00004 // Class:      HLTMuonDQMSource
00005 // 
00013 //
00014 // Original Author:  Muriel VANDER DONCKT *:0
00015 //         Created:  Wed Dec 12 09:55:42 CET 2007
00016 // $Id: HLTMuonDQMSource.cc,v 1.9 2008/10/16 16:35:15 hdyoo Exp $
00017 //
00018 //
00019 
00020 
00021 
00022 #include "DQM/HLTEvF/interface/HLTMuonDQMSource.h"
00023 #include "DQMServices/Core/interface/DQMStore.h"
00024 
00025 #include "DataFormats/Common/interface/Handle.h"
00026 #include "FWCore/Framework/interface/Event.h"
00027 #include "FWCore/Framework/interface/EventSetup.h"
00028 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00029 #include "FWCore/ServiceRegistry/interface/Service.h"
00030 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
00031 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
00032 #include "DataFormats/TrackReco/interface/Track.h"
00033 #include "DataFormats/MuonReco/interface/MuonFwd.h"
00034 #include "DataFormats/RecoCandidate/interface/IsoDeposit.h"
00035 #include "DataFormats/RecoCandidate/interface/IsoDepositFwd.h"
00036 #include "DataFormats/Common/interface/AssociationMap.h"
00037 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00038 #include "DataFormats/MuonSeed/interface/L2MuonTrajectorySeed.h"
00039 #include "DataFormats/MuonSeed/interface/L2MuonTrajectorySeedCollection.h"
00040 #include "DataFormats/TrajectorySeed/interface/TrajectorySeed.h"
00041 #include "DataFormats/L1Trigger/interface/L1MuonParticle.h"
00042 #include "DataFormats/L1Trigger/interface/L1MuonParticleFwd.h"
00043 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTCand.h"
00044 
00045 #include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeed.h"
00046 #include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeedCollection.h"
00047 
00048 #include "TMath.h" 
00049 
00050 
00051 using namespace std;
00052 using namespace edm;
00053 using namespace reco;
00054 using namespace l1extra;
00055 //
00056 // constructors and destructor
00057 //
00058 HLTMuonDQMSource::HLTMuonDQMSource( const edm::ParameterSet& ps ) :counterEvt_(0)
00059 
00060 {
00061   parameters_ = ps;
00062   verbose_ = parameters_.getUntrackedParameter < bool > ("verbose", false);
00063   monitorName_ = parameters_.getUntrackedParameter<string>("monitorName","HLT/HLTMonMuon");
00064   prescaleEvt_ = parameters_.getUntrackedParameter<int>("prescaleEvt", -1);
00065   coneSize_ = parameters_.getUntrackedParameter<double>("coneSize", 0.24);
00066   l2seedscollectionTag_ = parameters_.getUntrackedParameter<InputTag>("l2MuonSeedTag",edm::InputTag("hltL2MuonSeeds"));
00067   l2collectionTag_ = parameters_.getUntrackedParameter<InputTag>("l2MuonTag",edm::InputTag("hltL2MuonCandidates"));
00068   l3collectionTag_ = parameters_.getUntrackedParameter<InputTag>("l3MuonTag",edm::InputTag("hltL3MuonCandidates"));
00069   l2isolationTag_ = parameters_.getUntrackedParameter<InputTag>("l2IsolationTag",edm::InputTag("hltL2MuonIsolations"));
00070   l3isolationTag_ = parameters_.getUntrackedParameter<InputTag>("l3IsolationTag",edm::InputTag("hltL3MuonIsolations"));
00071 
00072    dbe_ = 0 ;
00073    dbe_ = Service < DQMStore > ().operator->();
00074    dbe_->setVerbose(0);
00075  
00076    outputFile_ =
00077        parameters_.getUntrackedParameter < std::string > ("outputFile", "");
00078    if (outputFile_.size() != 0) {
00079      LogWarning("HLTMuonDQMSource") << "Muon HLT Monitoring histograms will be saved to " 
00080                << outputFile_ << std::endl;
00081    }
00082    else {
00083      outputFile_ = "HLTMuonDQM.root";
00084    }
00085  
00086    bool disable =
00087      parameters_.getUntrackedParameter < bool > ("disableROOToutput", false);
00088    if (disable) {
00089      outputFile_ = "";
00090    }
00091 
00092    if (dbe_ != NULL) {
00093      dbe_->setCurrentFolder("HLT/HLTMonMuon");
00094    }
00095 
00096 
00097 }
00098 
00099 
00100 HLTMuonDQMSource::~HLTMuonDQMSource()
00101 {
00102    
00103   // do anything here that needs to be done at desctruction time
00104   // (e.g. close files, deallocate resources etc.)
00105   
00106 }
00107 
00108 
00109 //--------------------------------------------------------
00110 void HLTMuonDQMSource::beginJob(const EventSetup& context){
00111 
00112    if (dbe_) {
00113      dbe_->setCurrentFolder(monitorName_);
00114      dbe_->rmdir(monitorName_);
00115    }
00116  
00117  
00118    if (dbe_) {
00119      //dbe_->setCurrentFolder("monitorName_");
00120      if (monitorName_ != "" ) monitorName_ = monitorName_+"/" ;
00121      LogInfo("HLTMuonDQMSource") << "===>DQM event prescale = " << prescaleEvt_ << " events "<< endl;
00122      
00123      
00125      const int NBINS = 100; XMIN = 0; XMAX = 50;
00126 
00127      // create and cd into new folder
00128      char name[512], title[512];
00129      for ( int level = 1 ; level < 5 ; ++ level ) {
00130        sprintf(name,"Level%i",level);
00131        dbe_->setCurrentFolder(monitorName_+name);
00132        if (level==1)hl1quality = dbe_->book1D("h1L1Quality","GMT quality Flag", 8, 0., 8.);
00133        sprintf(name,"HLTMuonL%i_NMu",level);
00134        sprintf(title,"L%i number of muons",level);
00135        hNMu[level-1] = dbe_->book1D(name,title, 5, 0., 5.);
00136        hNMu[level-1]->setAxisTitle("Number of muons", 1);
00137        sprintf(name,"HLTMuonL%i_pt",level);
00138        sprintf(title,"L%i Pt",level);
00139        hpt[level-1] = dbe_->book1D(name,title, NBINS, 0., 100);
00140        hpt[level-1]->setAxisTitle("Pt", 1);
00141        sprintf(name,"HLTMuonL%i_eta",level);
00142        sprintf(title,"L%i Muon #eta",level);
00143        heta[level-1] = dbe_->book1D(name,title, NBINS, -2.5, 2.5);
00144        heta[level-1]->setAxisTitle("#eta", 1);
00145        sprintf(name,"HLTMuonL%i_phi",level);
00146        sprintf(title,"L%i Muon #phi",level);
00147        hphi[level-1] = dbe_->book1D(name,title, NBINS, -3.15, 3.15);
00148        hphi[level-1]->setAxisTitle("#phi", 1);
00149        sprintf(name,"HLTMuonL%i_etaphi",level);
00150        sprintf(title,"L%i Muon #eta vs #phi",level);
00151        hetaphi[level-1] = dbe_->book2D(name,title, NBINS, -3.15, 3.15,NBINS,-2.5, 2.5);
00152        hetaphi[level-1]->setAxisTitle("#phi", 1);
00153        hetaphi[level-1]->setAxisTitle("#eta", 2); 
00154        sprintf(name,"HLTMuonL%i_ptphi",level);
00155        sprintf(title,"L%i Muon pt vs #phi",level);         
00156        hptphi[level-1] = dbe_->book2D(name,title, NBINS, 0., 100.,NBINS,-3.15, 3.15);
00157        hptphi[level-1]->setAxisTitle("pt", 1);
00158        hptphi[level-1]->setAxisTitle("#phi", 2);
00159        sprintf(name,"HLTMuonL%i_pteta",level);
00160        sprintf(title,"L%i Muon pt vs #eta",level);         
00161        hpteta[level-1] = dbe_->book2D(name,title, NBINS, 0., 100.,NBINS,-2.5, 2.5);
00162        hpteta[level-1]->setAxisTitle("pt", 1);
00163        hpteta[level-1]->setAxisTitle("#eta", 2);
00164        sprintf(name,"HLTMuonL%i_nhit",level);
00165        sprintf(title,"L%i Number of Valid Hits",level);         
00166        hnhit[level-1] = dbe_->book1D(name,title, NBINS, 0., 100.);
00167        hnhit[level-1]->setAxisTitle("Number of Valid Hits", 1);
00168        sprintf(name,"HLTMuonL%i_charge",level);
00169        sprintf(title,"L%i Muon Charge",level);         
00170        hcharge[level-1]  = dbe_->book1D(name,title, 3, -1.5, 1.5);
00171        hcharge[level-1]->setAxisTitle("Charge", 1);
00172        if (level>1&&level<4){
00173          sprintf(name,"HLTMuonL%i_ptlx",level);
00174          sprintf(title,"L%i Muon 90 percent efficiency Pt",level);
00175          hptlx[level-2] = dbe_->book1D(name,title, NBINS, 0., 100);
00176          hptlx[level-2]->setAxisTitle("90% efficiency Pt", 1);
00177          sprintf(name,"HLTMuonL%i_dr",level);
00178          sprintf(title,"L%i Muon radial impact vs BeamSpot",level);         
00179          hdr[level-2] = dbe_->book1D(name,title, NBINS, -0.3, 0.3);
00180          hdr[level-2]->setAxisTitle("R Impact (cm) vs BeamSpot", 1);
00181          sprintf(name,"HLTMuonL%i_d0",level);
00182          sprintf(title,"L%i Muon radial impact vs (0,0)",level);         
00183          hd0[level-2] = dbe_->book1D(name,title, NBINS, -0.3, 0.3);
00184          hd0[level-2]->setAxisTitle("R Impact (cm) vs 0,0", 1);
00185          sprintf(name,"HLTMuonL%i_dz",level);
00186          sprintf(title,"L%i Muon Z impact",level);         
00187          hdz[level-2] = dbe_->book1D(name,title, NBINS, -25., 25.);
00188          hdz[level-2]->setAxisTitle("Z impact (cm)", 1);
00189          sprintf(name,"HLTMuonL%i_err0",level);
00190          sprintf(title,"L%i Muon Error on Pt",level);         
00191          herr0[level-2] = dbe_->book1D(name,title,NBINS, 0., 0.03);
00192          herr0[level-2]->setAxisTitle("Error on Pt", 1);
00193          sprintf(name,"HLTMuonL%i_iso",level);
00194          if (level==2)sprintf(title,"L%i Muon Energy in Isolation cone",level);         
00195          else if (level==3)sprintf(title,"L%i Muon SumPt in Isolation cone",level);               
00196          hiso[level-2]  = dbe_->book1D(name,title, NBINS, 0., 10./(level-2));
00197          if ( level==2)hiso[level-2]->setAxisTitle("Calo Energy in Iso Cone (GeV)", 1);
00198          else if ( level==3)hiso[level-2]->setAxisTitle("Sum Pt in Iso Cone (GeV)", 1);
00199          sprintf(name,"HLTMuonL%i_DiMuMass",level);
00200          sprintf(title,"L%i Opposite charge DiMuon invariant Mass",level);         
00201          hdimumass[level-2]= dbe_->book1D(name,title, NBINS, 0., 150.);
00202          hdimumass[level-2]->setAxisTitle("Di Muon Invariant Mass (GeV)");
00203          sprintf(name,"HLTMuonL%i_drphi",level);
00204          sprintf(title,"L%i #Deltar vs #phi",level);         
00205          hdrphi[level-2] = dbe_->bookProfile(name,title, NBINS, -3.15, 3.15,1,-999.,999.,"s");
00206          hdrphi[level-2]->setAxisTitle("#phi", 1);
00207          hdrphi[level-2]->setAxisTitle("#Deltar", 2);
00208          sprintf(name,"HLTMuonL%i_d0phi",level);
00209          sprintf(title,"L%i #Delta0 vs #phi",level);         
00210          hd0phi[level-2] = dbe_->bookProfile(name,title, NBINS, -3.15, 3.15,1,-999.,999.,"s");
00211          hd0phi[level-2]->setAxisTitle("#phi", 1);
00212          hd0phi[level-2]->setAxisTitle("#Delta0", 2);
00213          sprintf(name,"HLTMuonL%i_dzeta",level);
00214          sprintf(title,"L%i #Deltaz vs #eta",level);         
00215          hdzeta[level-2] = dbe_->bookProfile(name,title, NBINS,-2.5, 2.5,1,-999.,999.,"s");
00216          hdzeta[level-2]->setAxisTitle("#eta", 1);
00217          hdzeta[level-2]->setAxisTitle("#Deltaz", 2);
00218        }
00219        if (level < 3 ) {
00220          sprintf(name,"HLTMuonL%i_ptres",level);
00221          sprintf(title,"L%iMuon1/Pt - L%iMuon1/Pt",level,level+1);         
00222          hptres[level-1] = dbe_->book1D(name,title, NBINS, -0.1, 0.1);
00223          sprintf(title,"1/PtL%i - 1/PtL%i",level,level+1);         
00224          hptres[level-1]->setAxisTitle(title, 1);
00225          sprintf(name,"HLTMuonL%i_etares",level);
00226          sprintf(title,"L%i Muon #Delta#eta (wrt L%i)",level,level+1);         
00227          hetares[level-1] =dbe_->book1D(name,title, NBINS, -0.1, 0.1);
00228          hetares[level-1]->setAxisTitle("#Delta#eta", 1);
00229          sprintf(name,"HLTMuonL%i_phires",level);
00230          sprintf(title,"L%i Muon #Delta#phi (wrt L%i)",level,level+1);         
00231          hphires[level-1] =dbe_->book1D(name,title, NBINS, -0.1, 0.1);
00232          hphires[level-1]->setAxisTitle("#Delta#phi", 1);
00233          sprintf(name,"HLTMuonL%i_phiresphi",level);
00234          sprintf(title,"L%i Muon #Delta#phi vs #phi ",level);         
00235          hphiresphi[level-1] =dbe_->bookProfile(name,title, NBINS, -3.15, 3.15,1,-999.,999.,"s");
00236          hphiresphi[level-1]->setAxisTitle("<#Delta#phi>", 2);
00237          hphiresphi[level-1]->setAxisTitle("#phi", 1);
00238          sprintf(name,"HLTMuonL%i_etareseta",level);
00239          sprintf(title,"L%i Muon #Delta#eta vs #eta ",level);         
00240          hetareseta[level-1] =dbe_->bookProfile(name,title, NBINS,-2.5, 2.5,1,-999.,999.,"s");
00241          hetareseta[level-1]->setAxisTitle("<#Delta#eta>", 2);
00242          hetareseta[level-1]->setAxisTitle("#eta", 1);
00243          if (level  == 1 ){
00244            sprintf(name,"HLTMuonL%i3_ptres",level);
00245            sprintf(title,"L%iMuon1/Pt - L%iMuon1/Pt",level,level+2);         
00246            hptres[level+1] = dbe_->book1D(name,title, NBINS, -0.1, 0.1);
00247            sprintf(title,"1/PtL%i - 1/PtL%i",level,level+2);         
00248            hptres[level+1]->setAxisTitle(title, 1);
00249            sprintf(name,"HLTMuonL%i3_etares",level);
00250            sprintf(title,"L%i Muon #Delta#eta (wrt L%i)",level,level+2);         
00251            hetares[level+1] =dbe_->book1D(name,title, NBINS, -0.1, 0.1);
00252            hetares[level+1]->setAxisTitle("#Delta#eta", 1);
00253            sprintf(name,"HLTMuonL%i3_phires",level);
00254            sprintf(title,"L%i Muon #Delta#phi (wrt L%i)",level,level+2);         
00255            hphires[level+1] =dbe_->book1D(name,title, NBINS, -0.1, 0.1);
00256            hphires[level+1]->setAxisTitle("#Delta#phi", 1);
00257            sprintf(name,"HLTMuonL%i3_phiresphi",level);
00258            sprintf(title,"L%i Muon #Delta#phi vs #phi (wrt L3) ",level);         
00259            hphiresphi[level+1] =dbe_->bookProfile(name,title, NBINS, -3.15, 3.15,1,-999.,999.,"s");
00260            hphiresphi[level+1]->setAxisTitle("<#Delta#phi>", 2);
00261            hphiresphi[level+1]->setAxisTitle("#phi", 1);
00262            sprintf(name,"HLTMuonL%i3_etareseta",level);
00263            sprintf(title,"L%i Muon #Delta#eta vs #eta (wrt L3) ",level);         
00264            hetareseta[level+1] =dbe_->bookProfile(name,title, NBINS,-2.5, 2.5,1,-999.,999.,"s");
00265            hetareseta[level+1]->setAxisTitle("<#Delta#eta>", 2);
00266            hetareseta[level+1]->setAxisTitle("#eta", 1);
00267          }
00268        }
00269      }
00270      dbe_->showDirStructure();
00271   
00272      // Muon det id is 2 pushed in bits 28:31
00273      const unsigned int detector_id = 2<<28;
00274      dbe_->tagContents(monitorName_, detector_id);
00275    } 
00276 }
00277 
00278 //--------------------------------------------------------
00279 void HLTMuonDQMSource::beginRun(const edm::Run& r, const EventSetup& context) {
00280 
00281 }
00282 
00283 //--------------------------------------------------------
00284 void HLTMuonDQMSource::beginLuminosityBlock(const LuminosityBlock& lumiSeg, 
00285                                       const EventSetup& context) {
00286   
00287 }
00288 
00289 // ----------------------------------------------------------
00290 void HLTMuonDQMSource::analyze(const Event& iEvent, 
00291                          const EventSetup& iSetup )
00292 {  
00293   if ( !dbe_) return;
00294   counterEvt_++;
00295   if (prescaleEvt_ > 0 && counterEvt_%prescaleEvt_!=0) return;
00296   LogDebug("HLTMuonDQMSource") << " processing conterEvt_: " << counterEvt_ <<endl;
00297 
00298   Handle<RecoChargedCandidateCollection> l2mucands, l3mucands;
00299   iEvent.getByLabel (l2collectionTag_,l2mucands);
00300   iEvent.getByLabel (l3collectionTag_,l3mucands);
00301 
00302   RecoChargedCandidateCollection::const_iterator cand,cand2;
00303   Handle<L2MuonTrajectorySeedCollection> l2seeds; 
00304   iEvent.getByLabel (l2seedscollectionTag_,l2seeds);
00305 
00306 
00307   if (!l2seeds.failedToGet()) {
00308     hNMu[3]->Fill(l2seeds->size());
00309     L2MuonTrajectorySeedCollection::const_iterator l2seed;
00310     for (l2seed=l2seeds->begin() ; l2seed != l2seeds->end();++l2seed){
00311       PTrajectoryStateOnDet state=l2seed->startingState();
00312       float pt=state.parameters().momentum().perp();
00313       float eta=state.parameters().momentum().phi();
00314       float phi=state.parameters().momentum().eta();
00315       hcharge[3]->Fill(state.parameters().charge());
00316       hpt[3]->Fill(pt);
00317       hnhit[3]->Fill(l2seed->nHits());
00318       hphi[3]->Fill(phi);
00319       heta[3]->Fill(eta);
00320       hetaphi[3]->Fill(phi,eta);
00321       hptphi[3]->Fill(pt,phi);
00322       hpteta[3]->Fill(pt,eta);
00323       L1MuonParticleRef l1ref=l2seed->l1Particle();
00324       hcharge[0]->Fill(l1ref->charge());
00325       hpt[0]->Fill(l1ref->pt());
00326       hphi[0]->Fill(l1ref->phi());
00327       heta[0]->Fill(l1ref->eta());
00328       hetaphi[0]->Fill(l1ref->phi(),l1ref->eta());
00329       hptphi[0]->Fill(l1ref->pt(),l1ref->phi());
00330       hpteta[0]->Fill(l1ref->pt(),l1ref->eta());
00331       hl1quality->Fill(l1ref->gmtMuonCand().quality());
00332       if ( !l2mucands.failedToGet()) {
00333         for (cand=l2mucands->begin(); cand!=l2mucands->end(); ++cand) {
00334           TrackRef tk = cand->get<TrackRef>();
00335           RefToBase<TrajectorySeed> seed=tk->seedRef();
00336           if ( (l2seed->startingState()).detId() == (seed->startingState()).detId() ) {
00337             if(tk->pt()*l1ref->pt() != 0 )hptres[0]->Fill(1/tk->pt() - 1/l1ref->pt());
00338             hetares[0]->Fill(tk->eta()-l1ref->eta());
00339             hetareseta[0]->Fill(tk->eta(),tk->eta()-l1ref->eta());
00340             hphires[0]->Fill(tk->phi()-l1ref->phi());
00341             double dphi=tk->phi()-l1ref->phi();
00342             if (dphi>TMath::TwoPi())dphi-=2*TMath::TwoPi();
00343             else if (dphi<-TMath::TwoPi()) dphi+=TMath::TwoPi();
00344             hphiresphi[0]->Fill(tk->phi(),dphi);
00345             //find the L3 build from this L2
00346             
00347             if (!l3mucands.failedToGet()) {
00348               for (cand=l3mucands->begin(); cand!=l3mucands->end(); ++cand) {
00349                 TrackRef l3tk= cand->get<TrackRef>();
00350                 if (l3tk->seedRef().castTo<Ref<L3MuonTrajectorySeedCollection> >()->l2Track() == tk){
00351                   if(l1ref->pt()*l3tk->pt() != 0 )hptres[2]->Fill(1/l1ref->pt() - 1/l3tk->pt());
00352                   hetares[2]->Fill(l1ref->eta()-l3tk->eta());
00353                   hetareseta[2]->Fill(l1ref->eta(),l1ref->eta()-l3tk->eta());
00354                   hphires[2]->Fill(l1ref->phi()-l3tk->phi());
00355                   double dphi=l1ref->phi()-l3tk->phi();
00356                   if (dphi>TMath::TwoPi())dphi-=2*TMath::TwoPi();
00357                   else if (dphi<-TMath::TwoPi()) dphi+=TMath::TwoPi();
00358                   hphiresphi[2]->Fill(l3tk->phi(),dphi);
00359                   //break; //plot only once per L2?
00360                 }//if
00361               }//for
00362             }
00363             break;
00364           }
00365         }
00366       }
00367     }
00368   }
00369 
00370 
00371 
00372   reco::BeamSpot beamSpot;
00373   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00374   iEvent.getByLabel("offlineBeamSpot",recoBeamSpotHandle);
00375   if (!recoBeamSpotHandle.failedToGet())  beamSpot = *recoBeamSpotHandle;
00376 
00377   if (!l2mucands.failedToGet()) {
00378     LogDebug("HLTMuonDQMSource") << " filling L2 stuff " << endl;
00379     Handle<reco::IsoDepositMap> l2depMap;
00380     iEvent.getByLabel (l2isolationTag_,l2depMap);
00381     hNMu[1]->Fill(l2mucands->size());
00382     for (cand=l2mucands->begin(); cand!=l2mucands->end(); ++cand) {
00383       TrackRef tk = cand->get<TrackRef>();
00384       if (!l2depMap.failedToGet()) {
00385           LogDebug("HLTMuonDQMSource") << " filling L2 Iso stuff " << endl;
00386           if ( l2depMap->contains(tk.id()) ){
00387             reco::IsoDepositMap::value_type calDeposit= (*l2depMap)[tk];
00388             double dephlt = calDeposit.depositWithin(coneSize_);
00389             hiso[0]->Fill(dephlt);
00390           }
00391       }
00392     
00393       // eta cut
00394       hpt[1]->Fill(tk->pt());      
00395       double apar0 = fabs(tk->parameter(0));
00396       if (apar0>0)hptlx[0]->Fill((1+3.9*tk->error(0)/apar0)*tk->pt());      
00397       hcharge[1]->Fill(tk->charge()); 
00398       if ( tk->charge() != 0 ) {
00399         heta[1]->Fill(tk->eta());      
00400         hphi[1]->Fill(tk->phi()); 
00401         hetaphi[1]->Fill(tk->phi(),tk->eta()); 
00402         hptphi[1]->Fill(tk->pt(),tk->phi()); 
00403         hpteta[1]->Fill(tk->pt(),tk->eta()); 
00404         hnhit[1]->Fill(tk->numberOfValidHits()); 
00405         hd0[0]->Fill(tk->d0()); 
00406         if (!recoBeamSpotHandle.failedToGet()){
00407           hdr[0]->Fill(tk->dxy(beamSpot.position()));   
00408           hdrphi[0]->Fill(tk->phi(),tk->dxy(beamSpot.position())); 
00409         } 
00410         hd0phi[0]->Fill(tk->phi(),tk->d0()); 
00411         hdz[0]->Fill(tk->dz()); 
00412         hdzeta[0]->Fill(tk->eta(),tk->dz());
00413         herr0[0]->Fill(tk->error(0)); 
00414         cand2=cand;
00415         ++cand2;
00416         for (; cand2!=l2mucands->end(); cand2++) {
00417           TrackRef tk2=cand2->get<TrackRef>();
00418           if ( tk->charge()*tk2->charge() == -1 ){
00419             double mass=(cand->p4()+cand2->p4()).M();
00420             hdimumass[0]->Fill(mass);
00421           }
00422         }
00423       } else LogWarning("HLTMonMuon")<<"stop filling candidate with update@Vtx failure";
00424     }
00425   }
00426   if (!l3mucands.failedToGet()) {
00427     LogDebug("HLTMuonDQMSource") << " filling L3 stuff " << endl;
00428     hNMu[2]->Fill(l3mucands->size());
00429     Handle<reco::IsoDepositMap> l3depMap;
00430     iEvent.getByLabel (l3isolationTag_,l3depMap);
00431     for (cand=l3mucands->begin(); cand!=l3mucands->end(); ++cand) {
00432       TrackRef tk = cand->get<TrackRef>();
00433       if (!l3depMap.failedToGet()) {
00434         if ( l3depMap->contains(tk.id()) ){
00435           reco::IsoDepositMap::value_type calDeposit= (*l3depMap)[tk];
00436           double dephlt = calDeposit.depositWithin(coneSize_);
00437           hiso[1]->Fill(dephlt);
00438         }
00439       }
00440       // eta cut
00441       hpt[2]->Fill(tk->pt());      
00442       double apar0 = fabs(tk->parameter(0));
00443       if (apar0>0)hptlx[1]->Fill((1+2.2*tk->error(0)/apar0)*tk->pt());      
00444       heta[2]->Fill(tk->eta());      
00445       hphi[2]->Fill(tk->phi()); 
00446       hetaphi[2]->Fill(tk->phi(),tk->eta()); 
00447       hptphi[2]->Fill(tk->pt(),tk->phi()); 
00448       hpteta[2]->Fill(tk->pt(),tk->eta()); 
00449       hnhit[2]->Fill(tk->numberOfValidHits()); 
00450       hd0[1]->Fill(tk->d0()); 
00451       if (!recoBeamSpotHandle.failedToGet()) {
00452         hdr[1]->Fill(tk->dxy(beamSpot.position()));
00453         hdrphi[1]->Fill(tk->phi(),tk->dxy(beamSpot.position())); 
00454       }
00455       hd0phi[1]->Fill(tk->phi(),tk->d0()); 
00456       hdz[1]->Fill(tk->dz()); 
00457       hdzeta[1]->Fill(tk->eta(),tk->dz());
00458       herr0[1]->Fill(tk->error(0)); 
00459       hcharge[2]->Fill(tk->charge()); 
00460       cand2=cand;
00461       ++cand2;
00462 
00463       for (; cand2!=l3mucands->end(); cand2++) {
00464         TrackRef tk2=cand2->get<TrackRef>();
00465         if ( tk->charge()*tk2->charge() == -1 ){
00466           double mass=(cand->p4()+cand2->p4()).M();
00467           hdimumass[2]->Fill(mass);
00468         }
00469       }
00470       TrackRef l2tk = tk->seedRef().castTo<Ref<L3MuonTrajectorySeedCollection> >()->l2Track();
00471       if(tk->pt()*l2tk->pt() != 0 )hptres[1]->Fill(1/tk->pt() - 1/l2tk->pt());
00472       hetares[1]->Fill(tk->eta()-l2tk->eta());
00473       hetareseta[1]->Fill(tk->eta(),tk->eta()-l2tk->eta());
00474       hphires[1]->Fill(tk->phi()-l2tk->phi());
00475       double dphi=tk->phi()-l2tk->phi();
00476       if (dphi>TMath::TwoPi())dphi-=2*TMath::TwoPi();
00477       else if (dphi<-TMath::TwoPi()) dphi+=TMath::TwoPi();
00478       hphiresphi[1]->Fill(tk->phi(),dphi);
00479     }
00480   }  
00481 }
00482 
00483 
00484 
00485 
00486 //--------------------------------------------------------
00487 void HLTMuonDQMSource::endLuminosityBlock(const LuminosityBlock& lumiSeg, 
00488                                     const EventSetup& context) {
00489 }
00490 //--------------------------------------------------------
00491 void HLTMuonDQMSource::endRun(const Run& r, const EventSetup& context){
00492 }
00493 //--------------------------------------------------------
00494 void HLTMuonDQMSource::endJob(){
00495    LogInfo("HLTMonMuon") << "analyzed " << counterEvt_ << " events";
00496  
00497    //  if (outputFile_.size() != 0 && dbe_)
00498    //  dbe->save(outputFile_);
00499  
00500    return;
00501 }

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