CMS 3D CMS Logo

Classes | Public Member Functions | Protected Attributes

RecoMuonValidator Class Reference

#include <RecoMuonValidator.h>

Inheritance diagram for RecoMuonValidator:
edm::EDAnalyzer

List of all members.

Classes

struct  CommonME
struct  MuonME

Public Member Functions

virtual void analyze (const edm::Event &event, const edm::EventSetup &eventSetup)
virtual void beginRun (const edm::Run &, const edm::EventSetup &eventSetup)
virtual int countMuonHits (const reco::Track &track) const
virtual int countTrackerHits (const reco::Track &track) const
virtual void endRun ()
 RecoMuonValidator (const edm::ParameterSet &pset)
 ~RecoMuonValidator ()

Protected Attributes

CommonMEcommonME_
bool doAbsEta_
bool doAssoc_
TrackAssociatorBaseglbMuAssociator_
edm::InputTag glbMuAssocLabel_
edm::InputTag glbMuLabel_
MuonMEglbMuME_
edm::InputTag muonLabel_
std::string outputFileName_
edm::InputTag simLabel_
TrackAssociatorBasestaMuAssociator_
edm::InputTag staMuAssocLabel_
edm::InputTag staMuLabel_
MuonMEstaMuME_
std::string subDir_
DQMStoretheDQM
MuonServiceProxytheMuonService
TrackingParticleSelector tpSelector_
TrackAssociatorBasetrkMuAssociator_
edm::InputTag trkMuAssocLabel_
edm::InputTag trkMuLabel_
MuonMEtrkMuME_
unsigned int verbose_

Detailed Description

Definition at line 20 of file RecoMuonValidator.h.


Constructor & Destructor Documentation

RecoMuonValidator::RecoMuonValidator ( const edm::ParameterSet pset)

Definition at line 358 of file RecoMuonValidator.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), RecoMuonValidator::MuonME::hNSimHits_, edm::InputTag::instance(), edm::InputTag::label(), MuonServiceProxy_cff::MuonServiceProxy, HistoDimensions::nBinP, cmsCodeRules::cppFunctionSkipper::operator, and Pi.

{
  verbose_ = pset.getUntrackedParameter<unsigned int>("verbose", 0);

  outputFileName_ = pset.getUntrackedParameter<string>("outputFileName", "");

  // Set histogram dimensions
  HistoDimensions hDim;

  hDim.nBinP = pset.getUntrackedParameter<unsigned int>("nBinP");
  hDim.minP = pset.getUntrackedParameter<double>("minP");
  hDim.maxP = pset.getUntrackedParameter<double>("maxP");

  hDim.nBinPt = pset.getUntrackedParameter<unsigned int>("nBinPt");
  hDim.minPt = pset.getUntrackedParameter<double>("minPt");
  hDim.maxPt = pset.getUntrackedParameter<double>("maxPt");

  doAbsEta_ = pset.getUntrackedParameter<bool>("doAbsEta");
  hDim.doAbsEta = doAbsEta_;
  hDim.nBinEta  = pset.getUntrackedParameter<unsigned int>("nBinEta");
  hDim.minEta = pset.getUntrackedParameter<double>("minEta");
  hDim.maxEta = pset.getUntrackedParameter<double>("maxEta");

  hDim.nBinPhi  = pset.getUntrackedParameter<unsigned int>("nBinPhi");
  hDim.minPhi = pset.getUntrackedParameter<double>("minPhi", -TMath::Pi());
  hDim.maxPhi = pset.getUntrackedParameter<double>("maxPhi",  TMath::Pi());

  hDim.nBinErr  = pset.getUntrackedParameter<unsigned int>("nBinErr");
  hDim.nBinPull = pset.getUntrackedParameter<unsigned int>("nBinPull");

  hDim.wPull = pset.getUntrackedParameter<double>("wPull");

  hDim.minErrP = pset.getUntrackedParameter<double>("minErrP");
  hDim.maxErrP = pset.getUntrackedParameter<double>("maxErrP");

  hDim.minErrPt = pset.getUntrackedParameter<double>("minErrPt");
  hDim.maxErrPt = pset.getUntrackedParameter<double>("maxErrPt");

  hDim.minErrQPt = pset.getUntrackedParameter<double>("minErrQPt");
  hDim.maxErrQPt = pset.getUntrackedParameter<double>("maxErrQPt");

  hDim.minErrEta = pset.getUntrackedParameter<double>("minErrEta");
  hDim.maxErrEta = pset.getUntrackedParameter<double>("maxErrEta");

  hDim.minErrPhi = pset.getUntrackedParameter<double>("minErrPhi");
  hDim.maxErrPhi = pset.getUntrackedParameter<double>("maxErrPhi");

  hDim.minErrDxy = pset.getUntrackedParameter<double>("minErrDxy");
  hDim.maxErrDxy = pset.getUntrackedParameter<double>("maxErrDxy");

  hDim.minErrDz  = pset.getUntrackedParameter<double>("minErrDz" );
  hDim.maxErrDz  = pset.getUntrackedParameter<double>("maxErrDz" );

  hDim.nTrks = pset.getUntrackedParameter<unsigned int>("nTrks");
  hDim.nAssoc = pset.getUntrackedParameter<unsigned int>("nAssoc");
  hDim.nDof = pset.getUntrackedParameter<unsigned int>("nDof", 55);

  // Labels for simulation and reconstruction tracks
  simLabel_  = pset.getParameter<InputTag>("simLabel" );
  trkMuLabel_ = pset.getParameter<InputTag>("trkMuLabel");
  staMuLabel_ = pset.getParameter<InputTag>("staMuLabel");
  glbMuLabel_ = pset.getParameter<InputTag>("glbMuLabel");
  muonLabel_ = pset.getParameter<InputTag>("muonLabel");

  // Labels for sim-reco association
  doAssoc_ = pset.getUntrackedParameter<bool>("doAssoc", true);
  trkMuAssocLabel_ = pset.getParameter<InputTag>("trkMuAssocLabel");
  staMuAssocLabel_ = pset.getParameter<InputTag>("staMuAssocLabel");
  glbMuAssocLabel_ = pset.getParameter<InputTag>("glbMuAssocLabel");

//  seedPropagatorName_ = pset.getParameter<string>("SeedPropagator");

  ParameterSet tpset = pset.getParameter<ParameterSet>("tpSelector");
  tpSelector_ = TrackingParticleSelector(tpset.getParameter<double>("ptMin"),
                                         tpset.getParameter<double>("minRapidity"),
                                         tpset.getParameter<double>("maxRapidity"),
                                         tpset.getParameter<double>("tip"),
                                         tpset.getParameter<double>("lip"),
                                         tpset.getParameter<int>("minHit"),
                                         tpset.getParameter<bool>("signalOnly"),
                                         tpset.getParameter<bool>("chargedOnly"),
                                         tpset.getParameter<std::vector<int> >("pdgId"));

  // the service parameters
  ParameterSet serviceParameters 
    = pset.getParameter<ParameterSet>("ServiceParameters");
  theMuonService = new MuonServiceProxy(serviceParameters);

  // retrieve the instance of DQMService
  theDQM = 0;
  theDQM = Service<DQMStore>().operator->();

  if ( ! theDQM ) {
    LogError("RecoMuonValidator") << "DQMService not initialized\n";
    return;
  }

  subDir_ = pset.getUntrackedParameter<string>("subDir");
  if ( subDir_.empty() ) subDir_ = "RecoMuonV";
  if ( subDir_[subDir_.size()-1] == '/' ) subDir_.erase(subDir_.size()-1);

  // book histograms
  theDQM->cd();

  theDQM->setCurrentFolder(subDir_+"/Muons");

  commonME_ = new CommonME;
  trkMuME_ = new MuonME;
  staMuME_ = new MuonME;
  glbMuME_ = new MuonME;

  commonME_->hSimP_   = theDQM->book1D("SimP"  , "p of simTracks"    , hDim.nBinP  , hDim.minP  , hDim.maxP  );
  commonME_->hSimPt_  = theDQM->book1D("SimPt" , "p_{T} of simTracks", hDim.nBinPt , hDim.minPt , hDim.maxPt );
  commonME_->hSimEta_ = theDQM->book1D("SimEta", "#eta of simTracks" , hDim.nBinEta, hDim.minEta, hDim.maxEta);
  commonME_->hSimPhi_ = theDQM->book1D("SimPhi", "#phi of simTracks" , hDim.nBinPhi, hDim.minPhi, hDim.maxPhi);

  commonME_->hNSim_  = theDQM->book1D("NSim" , "Number of particles per event", hDim.nTrks, 0, hDim.nTrks);
  commonME_->hNMuon_ = theDQM->book1D("NMuon", "Number of muons per event"    , hDim.nTrks, 0, hDim.nTrks);

  const int nHits = 201;
  commonME_->hNSimHits_ = theDQM->book1D("NSimHits", "Number of simHits", nHits, -100.5, 100.5);

  commonME_->hTrkToGlbDiffNTrackerHits_ = theDQM->book1D("TrkGlbDiffNTrackerHits", "Difference of number of tracker hits (tkMuon - globalMuon)", nHits, -100.5, 100.5);
  commonME_->hStaToGlbDiffNMuonHits_ = theDQM->book1D("StaGlbDiffNMuonHits", "Difference of number of muon hits (staMuon - globalMuon)", nHits, -100.5, 100.5);

  commonME_->hTrkToGlbDiffNTrackerHitsEta_ = theDQM->book2D("TrkGlbDiffNTrackerHitsEta", "Difference of number of tracker hits (tkMuon - globalMuon)",   hDim.nBinEta, hDim.minEta, hDim.maxEta,nHits, -100.5, 100.5);
  commonME_->hStaToGlbDiffNMuonHitsEta_ = theDQM->book2D("StaGlbDiffNMuonHitsEta", "Difference of number of muon hits (staMuon - globalMuon)",   hDim.nBinEta, hDim.minEta, hDim.maxEta,nHits, -100.5, 100.5);

  commonME_->hTrkToGlbDiffNTrackerHitsPt_ = theDQM->book2D("TrkGlbDiffNTrackerHitsPt", "Difference of number of tracker hits (tkMuon - globalMuon)",  hDim.nBinPt, hDim.minPt, hDim.maxPt,nHits, -100.5, 100.5);
  commonME_->hStaToGlbDiffNMuonHitsPt_ = theDQM->book2D("StaGlbDiffNMuonHitsPt", "Difference of number of muon hits (staMuon - globalMuon)",  hDim.nBinPt, hDim.minPt, hDim.maxPt,nHits, -100.5, 100.5);

  // - histograms on tracking variables
  theDQM->setCurrentFolder(subDir_+"/Trk");
  theDQM->bookString("TrackLabel", trkMuLabel_.label()+"_"+trkMuLabel_.instance());
  theDQM->bookString("AssocLabel", trkMuAssocLabel_.label());

  theDQM->setCurrentFolder(subDir_+"/Sta");
  theDQM->bookString("TrackLabel", staMuLabel_.label()+"_"+staMuLabel_.instance());
  theDQM->bookString("AssocLabel", staMuAssocLabel_.label());

  theDQM->setCurrentFolder(subDir_+"/Glb");
  theDQM->bookString("TrackLabel", glbMuLabel_.label()+"_"+glbMuLabel_.instance());
  theDQM->bookString("AssocLabel", glbMuAssocLabel_.label());
  
  trkMuME_->bookHistograms(theDQM, subDir_+"/Trk", hDim);
  staMuME_->bookHistograms(theDQM, subDir_+"/Sta", hDim);
  glbMuME_->bookHistograms(theDQM, subDir_+"/Glb", hDim);

  if ( verbose_ > 0 ) theDQM->showDirStructure();

}
RecoMuonValidator::~RecoMuonValidator ( )

Definition at line 510 of file RecoMuonValidator.cc.

{
  if ( theMuonService ) delete theMuonService;
}

Member Function Documentation

void RecoMuonValidator::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
) [virtual]

Implements edm::EDAnalyzer.

Definition at line 542 of file RecoMuonValidator.cc.

References edm::View< T >::begin(), edm::AssociationMap< Tag >::end(), edm::View< T >::end(), reco::TrackBase::fill(), edm::AssociationMap< Tag >::find(), edm::Ref< C, T, F >::get(), i, edm::Handle< T >::product(), and edm::View< T >::size().

{
  if ( ! theDQM ) {
    LogError("RecoMuonValidator") << "DQMService not initialized\n";
    return;
  }

  // Get TrackingParticles
  Handle<TrackingParticleCollection> simHandle;
  event.getByLabel(simLabel_, simHandle);
  const TrackingParticleCollection simColl = *(simHandle.product());

  // Get Muon Tracks
  Handle<View<Track> > trkMuHandle;
  event.getByLabel(trkMuLabel_, trkMuHandle);
  View<Track> trkMuColl = *(trkMuHandle.product());

  Handle<View<Track> > staMuHandle;
  event.getByLabel(staMuLabel_, staMuHandle);
  View<Track> staMuColl = *(staMuHandle.product());

  Handle<View<Track> > glbMuHandle;
  event.getByLabel(glbMuLabel_, glbMuHandle);
  View<Track> glbMuColl = *(glbMuHandle.product());

  // Get Muons
  Handle<View<Muon> > muonHandle;
  event.getByLabel(muonLabel_, muonHandle);
  View<Muon> muonColl = *(muonHandle.product());

  // Get Association maps
  SimToRecoCollection simToTrkMuColl;
  SimToRecoCollection simToStaMuColl;
  SimToRecoCollection simToGlbMuColl;

  RecoToSimCollection trkMuToSimColl;
  RecoToSimCollection staMuToSimColl;
  RecoToSimCollection glbMuToSimColl;

  if ( doAssoc_ ) {
    // SimToReco associations
    simToTrkMuColl = trkMuAssociator_->associateSimToReco(trkMuHandle, simHandle, &event);
    simToStaMuColl = staMuAssociator_->associateSimToReco(staMuHandle, simHandle, &event);
    simToGlbMuColl = glbMuAssociator_->associateSimToReco(glbMuHandle, simHandle, &event);

    // // RecoToSim associations
    trkMuToSimColl = trkMuAssociator_->associateRecoToSim(trkMuHandle, simHandle, &event);
    staMuToSimColl = staMuAssociator_->associateRecoToSim(staMuHandle, simHandle, &event);
    glbMuToSimColl = glbMuAssociator_->associateRecoToSim(glbMuHandle, simHandle, &event);
  }
  else {
    // SimToReco associations
    Handle<SimToRecoCollection> simToTrkMuHandle;
    event.getByLabel(trkMuAssocLabel_, simToTrkMuHandle);
    simToTrkMuColl = *(simToTrkMuHandle.product());

    Handle<SimToRecoCollection> simToStaMuHandle;
    event.getByLabel(staMuAssocLabel_, simToStaMuHandle);
    simToStaMuColl = *(simToStaMuHandle.product());

    Handle<SimToRecoCollection> simToGlbMuHandle;
    event.getByLabel(glbMuAssocLabel_, simToGlbMuHandle);
    simToGlbMuColl = *(simToGlbMuHandle.product());

    // RecoToSim associations
    Handle<RecoToSimCollection> trkMuToSimHandle;
    event.getByLabel(trkMuAssocLabel_, trkMuToSimHandle);
    trkMuToSimColl = *(trkMuToSimHandle.product());

    Handle<RecoToSimCollection> staMuToSimHandle;
    event.getByLabel(staMuAssocLabel_, staMuToSimHandle);
    staMuToSimColl = *(staMuToSimHandle.product());

    Handle<RecoToSimCollection> glbMuToSimHandle;
    event.getByLabel(glbMuAssocLabel_, glbMuToSimHandle);
    glbMuToSimColl = *(glbMuToSimHandle.product());
  }

  const TrackingParticleCollection::size_type nSim = simColl.size();
  commonME_->hNSim_->Fill(nSim);

  commonME_->hNMuon_->Fill(muonColl.size());

  trkMuME_->hNTrks_->Fill(trkMuColl.size());
  staMuME_->hNTrks_->Fill(staMuColl.size());
  glbMuME_->hNTrks_->Fill(glbMuColl.size());

  // Analyzer reco::Muon
  for(View<Muon>::const_iterator iMuon = muonColl.begin();
      iMuon != muonColl.end(); ++iMuon) {
    int trkNTrackerHits = 0, glbNTrackerHits = 0;
    int staNMuonHits = 0, glbNMuonHits = 0;

    if ( iMuon->isTrackerMuon() ) {
      const TrackRef trkTrack = iMuon->track();

      trkNTrackerHits = countTrackerHits(*trkTrack);

      trkMuME_->hNTrackerHits_->Fill(trkNTrackerHits);
      trkMuME_->hNTrackerHits_vs_Pt_->Fill(trkTrack->pt(), trkNTrackerHits);
      trkMuME_->hNTrackerHits_vs_Eta_->Fill(trkTrack->eta(), trkNTrackerHits);
    }

    if ( iMuon->isStandAloneMuon() ) {
      const TrackRef staTrack = iMuon->standAloneMuon();

      staNMuonHits = countMuonHits(*staTrack);

      staMuME_->hNMuonHits_->Fill(staNMuonHits);
      staMuME_->hNMuonHits_vs_Pt_->Fill(staTrack->pt(), staNMuonHits);
      staMuME_->hNMuonHits_vs_Eta_->Fill(staTrack->eta(), staNMuonHits);

      staMuME_->hNTrksEta_->Fill(staTrack->eta());
      staMuME_->hNTrksPt_->Fill(staTrack->pt());
      
    }

    if ( iMuon->isGlobalMuon() ) {
      const TrackRef glbTrack = iMuon->combinedMuon();

      glbNTrackerHits = countTrackerHits(*glbTrack);
      glbNMuonHits = countMuonHits(*glbTrack);

      glbMuME_->hNTrackerHits_->Fill(glbNTrackerHits);
      glbMuME_->hNTrackerHits_vs_Pt_->Fill(glbTrack->pt(), glbNTrackerHits);
      glbMuME_->hNTrackerHits_vs_Eta_->Fill(glbTrack->eta(), glbNTrackerHits);

      glbMuME_->hNMuonHits_->Fill(glbNMuonHits);
      glbMuME_->hNMuonHits_vs_Pt_->Fill(glbTrack->pt(), glbNMuonHits);
      glbMuME_->hNMuonHits_vs_Eta_->Fill(glbTrack->eta(), glbNMuonHits);

      glbMuME_->hNTrksEta_->Fill(glbTrack->eta());
      glbMuME_->hNTrksPt_->Fill(glbTrack->pt());
      
      commonME_->hTrkToGlbDiffNTrackerHitsEta_->Fill(glbTrack->eta(),trkNTrackerHits-glbNTrackerHits);
      commonME_->hStaToGlbDiffNMuonHitsEta_->Fill(glbTrack->eta(),staNMuonHits-glbNMuonHits);
      
      commonME_->hTrkToGlbDiffNTrackerHitsPt_->Fill(glbTrack->pt(),trkNTrackerHits-glbNTrackerHits);
      commonME_->hStaToGlbDiffNMuonHitsPt_->Fill(glbTrack->pt(),staNMuonHits-glbNMuonHits);
      
      commonME_->hTrkToGlbDiffNTrackerHits_->Fill(trkNTrackerHits-glbNTrackerHits);
      commonME_->hStaToGlbDiffNMuonHits_->Fill(staNMuonHits-glbNMuonHits);

    }

  }

  // Analyzer reco::Track
  for(TrackingParticleCollection::size_type i=0; i<nSim; i++) {
    TrackingParticleRef simRef(simHandle, i);
    const TrackingParticle* simTP = simRef.get();
    if ( ! tpSelector_(*simTP) ) continue;

    const double simP   = simRef->p();
    const double simPt  = simRef->pt();
    const double simEta = doAbsEta_ ? fabs(simRef->eta()) : simRef->eta();
    const double simPhi = simRef->phi();

    const unsigned int nSimHits = simRef->pSimHit_end() - simRef->pSimHit_begin();

    commonME_->hSimP_  ->Fill(simP  );
    commonME_->hSimPt_ ->Fill(simPt );
    commonME_->hSimEta_->Fill(simEta);
    commonME_->hSimPhi_->Fill(simPhi);

    commonME_->hNSimHits_->Fill(nSimHits);

    // Get sim-reco association for a simRef
    vector<pair<RefToBase<Track>, double> > trkMuRefV, staMuRefV, glbMuRefV;
    if ( simToTrkMuColl.find(simRef) != simToTrkMuColl.end() ) {
      trkMuRefV = simToTrkMuColl[simRef];

      trkMuME_->hNSimToReco_->Fill(trkMuRefV.size());
      if ( ! trkMuRefV.empty() ) {
        const Track* trkMuTrack = trkMuRefV.begin()->first.get();
//        const double assocQuality = trkMuRefV.begin()->second;

        trkMuME_->fill(simTP, trkMuTrack);
      }
    }

    if ( simToStaMuColl.find(simRef) != simToStaMuColl.end() ) {
      staMuRefV = simToStaMuColl[simRef];

      staMuME_->hNSimToReco_->Fill(staMuRefV.size());
      if ( ! staMuRefV.empty() ) {
        const Track* staMuTrack = staMuRefV.begin()->first.get();
//        const double assocQuality = staMuRefV.begin().second;

        staMuME_->fill(simTP, staMuTrack);
      }
    }

    if ( simToGlbMuColl.find(simRef) != simToGlbMuColl.end() ) {
      glbMuRefV = simToGlbMuColl[simRef];

      glbMuME_->hNSimToReco_->Fill(glbMuRefV.size());
      if ( ! glbMuRefV.empty() ) {
        const Track* glbMuTrack = glbMuRefV.begin()->first.get();
//        const double assocQuality = glbMuRefV.begin().second;

        glbMuME_->fill(simTP, glbMuTrack);
      }
    }
  }
}
void RecoMuonValidator::beginRun ( const edm::Run ,
const edm::EventSetup eventSetup 
) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 515 of file RecoMuonValidator.cc.

References edm::EventSetup::get(), and edm::ESHandle< T >::product().

{
  if ( theMuonService ) theMuonService->update(eventSetup);

  trkMuAssociator_ = 0;
  staMuAssociator_ = 0;
  glbMuAssociator_ = 0;
  if ( doAssoc_ ) {
    ESHandle<TrackAssociatorBase> trkMuAssocHandle;
    eventSetup.get<TrackAssociatorRecord>().get(trkMuAssocLabel_.label(), trkMuAssocHandle);
    trkMuAssociator_ = const_cast<TrackAssociatorBase*>(trkMuAssocHandle.product());

    ESHandle<TrackAssociatorBase> staMuAssocHandle;
    eventSetup.get<TrackAssociatorRecord>().get(staMuAssocLabel_.label(), staMuAssocHandle);
    staMuAssociator_ = const_cast<TrackAssociatorBase*>(staMuAssocHandle.product());

    ESHandle<TrackAssociatorBase> glbMuAssocHandle;
    eventSetup.get<TrackAssociatorRecord>().get(glbMuAssocLabel_.label(), glbMuAssocHandle);
    glbMuAssociator_ = const_cast<TrackAssociatorBase*>(glbMuAssocHandle.product());
  }
}
int RecoMuonValidator::countMuonHits ( const reco::Track track) const [virtual]

Definition at line 750 of file RecoMuonValidator.cc.

References prof2calltree::count, DetId::det(), DetId::Muon, reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), and query::result.

                                                             {
  TransientTrackingRecHit::ConstRecHitContainer result;
  
  int count = 0;

  for (trackingRecHit_iterator hit = track.recHitsBegin(); hit != track.recHitsEnd(); ++hit) {
    if((*hit)->isValid()) {
      DetId recoid = (*hit)->geographicalId();
      if ( recoid.det() == DetId::Muon ) count++;
    }
  }
  return count;
}
int RecoMuonValidator::countTrackerHits ( const reco::Track track) const [virtual]

Definition at line 765 of file RecoMuonValidator.cc.

References prof2calltree::count, DetId::det(), reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), query::result, and align::Tracker.

                                                                {
  TransientTrackingRecHit::ConstRecHitContainer result;
  
  int count = 0;

  for (trackingRecHit_iterator hit = track.recHitsBegin(); hit != track.recHitsEnd(); ++hit) {
    if((*hit)->isValid()) {
      DetId recoid = (*hit)->geographicalId();
      if ( recoid.det() == DetId::Tracker ) count++;
    }
  }
  return count;
}
void RecoMuonValidator::endRun ( void  ) [virtual]

Definition at line 537 of file RecoMuonValidator.cc.

{
  if ( theDQM && ! outputFileName_.empty() ) theDQM->save(outputFileName_);
}

Member Data Documentation

Definition at line 60 of file RecoMuonValidator.h.

bool RecoMuonValidator::doAbsEta_ [protected]

Definition at line 51 of file RecoMuonValidator.h.

bool RecoMuonValidator::doAssoc_ [protected]

Definition at line 52 of file RecoMuonValidator.h.

Definition at line 55 of file RecoMuonValidator.h.

Definition at line 43 of file RecoMuonValidator.h.

Definition at line 38 of file RecoMuonValidator.h.

Definition at line 57 of file RecoMuonValidator.h.

Definition at line 39 of file RecoMuonValidator.h.

std::string RecoMuonValidator::outputFileName_ [protected]

Definition at line 45 of file RecoMuonValidator.h.

Definition at line 35 of file RecoMuonValidator.h.

Definition at line 55 of file RecoMuonValidator.h.

Definition at line 42 of file RecoMuonValidator.h.

Definition at line 37 of file RecoMuonValidator.h.

Definition at line 57 of file RecoMuonValidator.h.

std::string RecoMuonValidator::subDir_ [protected]

Definition at line 46 of file RecoMuonValidator.h.

Definition at line 49 of file RecoMuonValidator.h.

Definition at line 48 of file RecoMuonValidator.h.

Definition at line 54 of file RecoMuonValidator.h.

Definition at line 55 of file RecoMuonValidator.h.

Definition at line 41 of file RecoMuonValidator.h.

Definition at line 36 of file RecoMuonValidator.h.

Definition at line 57 of file RecoMuonValidator.h.

unsigned int RecoMuonValidator::verbose_ [protected]

Definition at line 33 of file RecoMuonValidator.h.