#include <DQM/TrackingMonitor/src/TrackSplittingMonitor.cc>
Monitoring source for general quantities related to tracks.
Definition at line 36 of file TrackSplittingMonitor.h.
TrackSplittingMonitor::TrackSplittingMonitor | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 37 of file TrackSplittingMonitor.cc.
References conf_, dqmStore_, and cppFunctionSkipper::operator.
{ dqmStore_ = edm::Service<DQMStore>().operator->(); conf_ = iConfig; }
TrackSplittingMonitor::~TrackSplittingMonitor | ( | ) |
Definition at line 42 of file TrackSplittingMonitor.cc.
{ }
void TrackSplittingMonitor::analyze | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 159 of file TrackSplittingMonitor.cc.
References cscGeometry, reco::TrackBase::d0(), d0Cut_, reco::TrackBase::d0Error(), ddxyAbsoluteResiduals_global_, ddxyAbsoluteResiduals_tracker_, ddxyNormalizedResiduals_global_, ddxyNormalizedResiduals_tracker_, dtGeometry, reco::TrackBase::dz(), dzCut_, reco::TrackBase::dzError(), MonitorElement::Fill(), edm::EventSetup::get(), edm::Event::getByLabel(), reco::Muon::globalTrack(), edm::HandleBase::isValid(), norchiCut_, reco::TrackBase::normalizedChi2(), reco::TrackBase::phi(), reco::TrackBase::phiError(), PixelSubdetector::PixelBarrel, pixelHitsPerLeg_, plotMuons_, funct::pow(), reco::TrackBase::pt(), ptCut_, reco::TrackBase::ptError(), reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), rpcGeometry, RecoMuonCosmics_cff::splitMuons, splitMuons_, splitTracks_, mathSSE::sqrt(), theGeometry, theMagField, reco::TrackBase::theta(), reco::TrackBase::thetaError(), and totalHitsPerLeg_.
{ iSetup.get<IdealMagneticFieldRecord>().get(theMagField); iSetup.get<TrackerDigiGeometryRecord>().get(theGeometry); iSetup.get<MuonGeometryRecord>().get(dtGeometry); iSetup.get<MuonGeometryRecord>().get(cscGeometry); iSetup.get<MuonGeometryRecord>().get(rpcGeometry); edm::Handle<std::vector<reco::Track> > splitTracks; iEvent.getByLabel(splitTracks_, splitTracks); if (!splitTracks.isValid()) return; edm::Handle<std::vector<reco::Muon> > splitMuons; if (plotMuons_){ iEvent.getByLabel(splitMuons_, splitMuons); } if (splitTracks->size() == 2){ // check that there are 2 tracks in split track collection edm::LogInfo("TrackSplittingMonitor") << "Split Track size: " << splitTracks->size(); // split tracks calculations reco::Track track1 = splitTracks->at(0); reco::Track track2 = splitTracks->at(1); // -------------------------- basic selection --------------------------- // hit counting // looping through the hits for track 1 double nRechits1 =0; double nRechitinBPIX1 =0; for (trackingRecHit_iterator iHit = track1.recHitsBegin(); iHit != track1.recHitsEnd(); ++iHit) { if((*iHit)->isValid()) { nRechits1++; int type =(*iHit)->geographicalId().subdetId(); if(type==int(PixelSubdetector::PixelBarrel)){++nRechitinBPIX1;} } } // looping through the hits for track 2 double nRechits2 =0; double nRechitinBPIX2 =0; for (trackingRecHit_iterator iHit = track2.recHitsBegin(); iHit != track2.recHitsEnd(); ++iHit) { if((*iHit)->isValid()) { nRechits2++; int type =(*iHit)->geographicalId().subdetId(); if(type==int(PixelSubdetector::PixelBarrel)){++nRechitinBPIX2;} } } // DCA of each track double d01 = track1.d0(); double dz1 = track1.dz(); double d02 = track2.d0(); double dz2 = track2.dz(); // pT of each track double pt1 = track1.pt(); double pt2 = track2.pt(); // chi2 of each track double norchi1 = track1.normalizedChi2(); double norchi2 = track2.normalizedChi2(); // basic selection // pixel hits and total hits if ((nRechitinBPIX1 >= pixelHitsPerLeg_)&&(nRechitinBPIX1 >= pixelHitsPerLeg_)&&(nRechits1 >= totalHitsPerLeg_)&&(nRechits2 >= totalHitsPerLeg_)){ // dca cut if ( ((fabs(d01) < d0Cut_))&&(fabs(d02) < d0Cut_)&&(fabs(dz2) < dzCut_)&&(fabs(dz2) < dzCut_) ){ // pt cut if ( (pt1+pt2)/2 < ptCut_){ // chi2 cut if ( (norchi1 < norchiCut_)&&(norchi2 < norchiCut_) ){ // passed all cuts... edm::LogInfo("TrackSplittingMonitor") << " Setected after all cuts ?"; double ddxyVal = d01 - d02; double ddzVal = dz1 - dz2; double dphiVal = track1.phi() - track2.phi(); double dthetaVal = track1.theta() - track2.theta(); double dptVal = pt1 - pt2; double dcurvVal = (1/pt1) - (1/pt2); double d01ErrVal = track1.d0Error(); double d02ErrVal = track2.d0Error(); double dz1ErrVal = track1.dzError(); double dz2ErrVal = track2.dzError(); double phi1ErrVal = track1.phiError(); double phi2ErrVal = track2.phiError(); double theta1ErrVal = track1.thetaError(); double theta2ErrVal = track2.thetaError(); double pt1ErrVal = track1.ptError(); double pt2ErrVal = track2.ptError(); ddxyAbsoluteResiduals_tracker_->Fill( 10000.0*ddxyVal/sqrt(2.0) ); ddxyAbsoluteResiduals_tracker_->Fill( 10000.0*ddzVal/sqrt(2.0) ); ddxyAbsoluteResiduals_tracker_->Fill( 1000.0*dphiVal/sqrt(2.0) ); ddxyAbsoluteResiduals_tracker_->Fill( 1000.0*dthetaVal/sqrt(2.0) ); ddxyAbsoluteResiduals_tracker_->Fill( dptVal/sqrt(2.0) ); ddxyAbsoluteResiduals_tracker_->Fill( dcurvVal/sqrt(2.0) ); ddxyNormalizedResiduals_tracker_->Fill( ddxyVal/sqrt( d01ErrVal*d01ErrVal + d02ErrVal*d02ErrVal ) ); ddxyNormalizedResiduals_tracker_->Fill( ddzVal/sqrt( dz1ErrVal*dz1ErrVal + dz2ErrVal*dz2ErrVal ) ); ddxyNormalizedResiduals_tracker_->Fill( dphiVal/sqrt( phi1ErrVal*phi1ErrVal + phi2ErrVal*phi2ErrVal ) ); ddxyNormalizedResiduals_tracker_->Fill( dthetaVal/sqrt( theta1ErrVal*theta1ErrVal + theta2ErrVal*theta2ErrVal ) ); ddxyNormalizedResiduals_tracker_->Fill( dptVal/sqrt( pt1ErrVal*pt1ErrVal + pt2ErrVal*pt2ErrVal ) ); ddxyNormalizedResiduals_tracker_->Fill( dcurvVal/sqrt( pow(pt1ErrVal,2)/pow(pt1,4) + pow(pt2ErrVal,2)/pow(pt2,4) ) ); // if do the same for split muons if (plotMuons_ && splitMuons.isValid()){ int gmCtr = 0; bool topGlobalMuonFlag = false; bool bottomGlobalMuonFlag = false; int topGlobalMuon = -1; int bottomGlobalMuon = -1; double topGlobalMuonNorchi2 = 1e10; double bottomGlobalMuonNorchi2 = 1e10; // check if usable split global muons for (std::vector<reco::Muon>::const_iterator gmI = splitMuons->begin(); gmI != splitMuons->end(); gmI++){ if ( gmI->isTrackerMuon() && gmI->isStandAloneMuon() && gmI->isGlobalMuon() ){ reco::TrackRef trackerTrackRef1( splitTracks, 0 ); reco::TrackRef trackerTrackRef2( splitTracks, 1 ); if (gmI->innerTrack() == trackerTrackRef1){ if (gmI->globalTrack()->normalizedChi2() < topGlobalMuonNorchi2){ topGlobalMuonFlag = true; topGlobalMuonNorchi2 = gmI->globalTrack()->normalizedChi2(); topGlobalMuon = gmCtr; } } if (gmI->innerTrack() == trackerTrackRef2){ if (gmI->globalTrack()->normalizedChi2() < bottomGlobalMuonNorchi2){ bottomGlobalMuonFlag = true; bottomGlobalMuonNorchi2 = gmI->globalTrack()->normalizedChi2(); bottomGlobalMuon = gmCtr; } } } gmCtr++; } if (bottomGlobalMuonFlag && topGlobalMuonFlag) { reco::Muon muonTop = splitMuons->at( topGlobalMuon ); reco::Muon muonBottom = splitMuons->at( bottomGlobalMuon ); reco::TrackRef glb1 = muonTop.globalTrack(); reco::TrackRef glb2 = muonBottom.globalTrack(); double ddxyValGlb = glb1->d0() - glb2->d0(); double ddzValGlb = glb1->dz() - glb2->dz(); double dphiValGlb = glb1->phi() - glb2->phi(); double dthetaValGlb = glb1->theta() - glb2->theta(); double dptValGlb = glb1->pt() - glb2->pt(); double dcurvValGlb = (1/glb1->pt()) - (1/glb2->pt()); double d01ErrValGlb = glb1->d0Error(); double d02ErrValGlb = glb2->d0Error(); double dz1ErrValGlb = glb1->dzError(); double dz2ErrValGlb = glb2->dzError(); double phi1ErrValGlb = glb1->phiError(); double phi2ErrValGlb = glb2->phiError(); double theta1ErrValGlb = glb1->thetaError(); double theta2ErrValGlb = glb2->thetaError(); double pt1ErrValGlb = glb1->ptError(); double pt2ErrValGlb = glb2->ptError(); ddxyAbsoluteResiduals_global_->Fill( 10000.0*ddxyValGlb/sqrt(2.0) ); ddxyAbsoluteResiduals_global_->Fill( 10000.0*ddzValGlb/sqrt(2.0) ); ddxyAbsoluteResiduals_global_->Fill( 1000.0*dphiValGlb/sqrt(2.0) ); ddxyAbsoluteResiduals_global_->Fill( 1000.0*dthetaValGlb/sqrt(2.0) ); ddxyAbsoluteResiduals_global_->Fill( dptValGlb/sqrt(2.0) ); ddxyAbsoluteResiduals_global_->Fill( dcurvValGlb/sqrt(2.0) ); ddxyNormalizedResiduals_global_->Fill( ddxyValGlb/sqrt( d01ErrValGlb*d01ErrValGlb + d02ErrValGlb*d02ErrValGlb ) ); ddxyNormalizedResiduals_global_->Fill( ddzValGlb/sqrt( dz1ErrValGlb*dz1ErrValGlb + dz2ErrValGlb*dz2ErrValGlb ) ); ddxyNormalizedResiduals_global_->Fill( dphiValGlb/sqrt( phi1ErrValGlb*phi1ErrValGlb + phi2ErrValGlb*phi2ErrValGlb ) ); ddxyNormalizedResiduals_global_->Fill( dthetaValGlb/sqrt( theta1ErrValGlb*theta1ErrValGlb + theta2ErrValGlb*theta2ErrValGlb ) ); ddxyNormalizedResiduals_global_->Fill( dptValGlb/sqrt( pt1ErrValGlb*pt1ErrValGlb + pt2ErrValGlb*pt2ErrValGlb ) ); ddxyNormalizedResiduals_global_->Fill( dcurvValGlb/sqrt( pow(pt1ErrValGlb,2)/pow(pt1,4) + pow(pt2ErrValGlb,2)/pow(pt2,4) ) ); } } // end of split muons loop } } } } } }
void TrackSplittingMonitor::beginJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 45 of file TrackSplittingMonitor.cc.
References DQMStore::book1D(), conf_, d0Cut_, dcurvAbsoluteResiduals_global_, dcurvAbsoluteResiduals_tracker_, dcurvNormalizedResiduals_global_, dcurvNormalizedResiduals_tracker_, ddxyAbsoluteResiduals_global_, ddxyAbsoluteResiduals_tracker_, ddxyNormalizedResiduals_global_, ddxyNormalizedResiduals_tracker_, ddzAbsoluteResiduals_global_, ddzAbsoluteResiduals_tracker_, ddzNormalizedResiduals_global_, ddzNormalizedResiduals_tracker_, dphiAbsoluteResiduals_global_, dphiAbsoluteResiduals_tracker_, dphiNormalizedResiduals_global_, dphiNormalizedResiduals_tracker_, dptAbsoluteResiduals_global_, dptAbsoluteResiduals_tracker_, dptNormalizedResiduals_global_, dptNormalizedResiduals_tracker_, dqmStore_, dthetaAbsoluteResiduals_global_, dthetaAbsoluteResiduals_tracker_, dthetaNormalizedResiduals_global_, dthetaNormalizedResiduals_tracker_, dzCut_, edm::ParameterSet::getParameter(), norchiCut_, pixelHitsPerLeg_, plotMuons_, ptCut_, MonitorElement::setAxisTitle(), DQMStore::setCurrentFolder(), splitMuons_, splitTracks_, and totalHitsPerLeg_.
{ std::string MEFolderName = conf_.getParameter<std::string>("FolderName"); dqmStore_->setCurrentFolder(MEFolderName); //get input tags splitTracks_ = conf_.getParameter< edm::InputTag >("splitTrackCollection"); splitMuons_ = conf_.getParameter< edm::InputTag >("splitMuonCollection"); plotMuons_ = conf_.getParameter<bool>("ifPlotMuons"); // cuts pixelHitsPerLeg_ = conf_.getParameter<int>("pixelHitsPerLeg"); totalHitsPerLeg_ = conf_.getParameter<int>("totalHitsPerLeg"); d0Cut_ = conf_.getParameter<double>("d0Cut"); dzCut_ = conf_.getParameter<double>("dzCut"); ptCut_ = conf_.getParameter<double>("ptCut"); norchiCut_ = conf_.getParameter<double>("norchiCut"); // bin declarations int ddxyBin = conf_.getParameter<int>("ddxyBin"); double ddxyMin = conf_.getParameter<double>("ddxyMin"); double ddxyMax = conf_.getParameter<double>("ddxyMax"); int ddzBin = conf_.getParameter<int>("ddzBin"); double ddzMin = conf_.getParameter<double>("ddzMin"); double ddzMax = conf_.getParameter<double>("ddzMax"); int dphiBin = conf_.getParameter<int>("dphiBin"); double dphiMin = conf_.getParameter<double>("dphiMin"); double dphiMax = conf_.getParameter<double>("dphiMax"); int dthetaBin = conf_.getParameter<int>("dthetaBin"); double dthetaMin = conf_.getParameter<double>("dthetaMin"); double dthetaMax = conf_.getParameter<double>("dthetaMax"); int dptBin = conf_.getParameter<int>("dptBin"); double dptMin = conf_.getParameter<double>("dptMin"); double dptMax = conf_.getParameter<double>("dptMax"); int dcurvBin = conf_.getParameter<int>("dcurvBin"); double dcurvMin = conf_.getParameter<double>("dcurvMin"); double dcurvMax = conf_.getParameter<double>("dcurvMax"); int normBin = conf_.getParameter<int>("normBin"); double normMin = conf_.getParameter<double>("normMin"); double normMax = conf_.getParameter<double>("normMax"); // declare histogram ddxyAbsoluteResiduals_tracker_ = dqmStore_->book1D( "ddxyAbsoluteResiduals_tracker", "ddxyAbsoluteResiduals_tracker", ddxyBin, ddxyMin, ddxyMax ); ddzAbsoluteResiduals_tracker_ = dqmStore_->book1D( "ddzAbsoluteResiduals_tracker", "ddzAbsoluteResiduals_tracker", ddzBin, ddzMin, ddzMax ); dphiAbsoluteResiduals_tracker_ = dqmStore_->book1D( "dphiAbsoluteResiduals_tracker", "dphiAbsoluteResiduals_tracker", dphiBin, dphiMin, dphiMax ); dthetaAbsoluteResiduals_tracker_ = dqmStore_->book1D( "dthetaAbsoluteResiduals_tracker", "dthetaAbsoluteResiduals_tracker", dthetaBin, dthetaMin, dthetaMax ); dptAbsoluteResiduals_tracker_ = dqmStore_->book1D( "dptAbsoluteResiduals_tracker", "dptAbsoluteResiduals_tracker", dptBin, dptMin, dptMax ); dcurvAbsoluteResiduals_tracker_ = dqmStore_->book1D( "dcurvAbsoluteResiduals_tracker", "dcurvAbsoluteResiduals_tracker", dcurvBin, dcurvMin, dcurvMax ); ddxyNormalizedResiduals_tracker_ = dqmStore_->book1D( "ddxyNormalizedResiduals_tracker", "ddxyNormalizedResiduals_tracker", normBin, normMin, normMax ); ddzNormalizedResiduals_tracker_ = dqmStore_->book1D( "ddzNormalizedResiduals_tracker", "ddzNormalizedResiduals_tracker", normBin, normMin, normMax ); dphiNormalizedResiduals_tracker_ = dqmStore_->book1D( "dphiNormalizedResiduals_tracker", "dphiNormalizedResiduals_tracker", normBin, normMin, normMax ); dthetaNormalizedResiduals_tracker_ = dqmStore_->book1D( "dthetaNormalizedResiduals_tracker", "dthetaNormalizedResiduals_tracker", normBin, normMin, normMax ); dptNormalizedResiduals_tracker_ = dqmStore_->book1D( "dptNormalizedResiduals_tracker", "dptNormalizedResiduals_tracker", normBin, normMin, normMax ); dcurvNormalizedResiduals_tracker_ = dqmStore_->book1D( "dcurvNormalizedResiduals_tracker", "dcurvNormalizedResiduals_tracker", normBin, normMin, normMax ); if (plotMuons_){ ddxyAbsoluteResiduals_global_ = dqmStore_->book1D( "ddxyAbsoluteResiduals_global", "ddxyAbsoluteResiduals_global", ddxyBin, ddxyMin, ddxyMax ); ddzAbsoluteResiduals_global_ = dqmStore_->book1D( "ddzAbsoluteResiduals_global", "ddzAbsoluteResiduals_global", ddzBin, ddzMin, ddzMax ); dphiAbsoluteResiduals_global_ = dqmStore_->book1D( "dphiAbsoluteResiduals_global", "dphiAbsoluteResiduals_global", dphiBin, dphiMin, dphiMax ); dthetaAbsoluteResiduals_global_ = dqmStore_->book1D( "dthetaAbsoluteResiduals_global", "dthetaAbsoluteResiduals_global", dthetaBin, dthetaMin, dthetaMax ); dptAbsoluteResiduals_global_ = dqmStore_->book1D( "dptAbsoluteResiduals_global", "dptAbsoluteResiduals_global", dptBin, dptMin, dptMax ); dcurvAbsoluteResiduals_global_ = dqmStore_->book1D( "dcurvAbsoluteResiduals_global", "dcurvAbsoluteResiduals_global", dcurvBin, dcurvMin, dcurvMax ); ddxyNormalizedResiduals_global_ = dqmStore_->book1D( "ddxyNormalizedResiduals_global", "ddxyNormalizedResiduals_global", normBin, normMin, normMax ); ddzNormalizedResiduals_global_ = dqmStore_->book1D( "ddzNormalizedResiduals_global", "ddzNormalizedResiduals_global", normBin, normMin, normMax ); dphiNormalizedResiduals_global_ = dqmStore_->book1D( "dphiNormalizedResiduals_global", "dphiNormalizedResiduals_global", normBin, normMin, normMax ); dthetaNormalizedResiduals_global_ = dqmStore_->book1D( "dthetaNormalizedResiduals_global", "dthetaNormalizedResiduals_global", normBin, normMin, normMax ); dptNormalizedResiduals_global_ = dqmStore_->book1D( "dptNormalizedResiduals_global", "dptNormalizedResiduals_global", normBin, normMin, normMax ); dcurvNormalizedResiduals_global_ = dqmStore_->book1D( "dcurvNormalizedResiduals_global", "dcurvNormalizedResiduals_global", normBin, normMin, normMax ); } ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta d_{xy})/#sqrt{2} [#mum]" ); ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta d_{z})/#sqrt{2} [#mum]" ); ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta #phi)/#sqrt{2} [mrad]" ); ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta #theta)/#sqrt{2} [mrad]" ); ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta pT)/#sqrt{2} [GeV]" ); ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta (1/pT))/#sqrt{2} [GeV^{-1}]" ); ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta d_{xy}/#sigma(d_{xy}" ); ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta d_{z}/#sigma(d_{z})" ); ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta #phi/#sigma(d_{#phi})" ); ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta #theta/#sigma(d_{#theta})" ); ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta p_{T}/#sigma(p_{T})" ); ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta 1/p_{T}/#sigma(1/p_{T})" ); if (plotMuons_){ ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta d_{xy})/#sqrt{2} [#mum]" ); ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta d_{z})/#sqrt{2} [#mum]" ); ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta #phi)/#sqrt{2} [mrad]" ); ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta #theta)/#sqrt{2} [mrad]" ); ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta pT)/#sqrt{2} [GeV]" ); ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta (1/pT))/#sqrt{2} [GeV^{-1}]" ); ddxyNormalizedResiduals_global_->setAxisTitle( "#delta d_{xy}/#sigma(d_{xy}" ); ddxyNormalizedResiduals_global_->setAxisTitle( "#delta d_{z}/#sigma(d_{z})" ); ddxyNormalizedResiduals_global_->setAxisTitle( "#delta #phi/#sigma(d_{#phi})" ); ddxyNormalizedResiduals_global_->setAxisTitle( "#delta #theta/#sigma(d_{#theta})" ); ddxyNormalizedResiduals_global_->setAxisTitle( "#delta p_{T}/#sigma(p_{T})" ); ddxyNormalizedResiduals_global_->setAxisTitle( "#delta 1/p_{T}/#sigma(1/p_{T})" ); } }
void TrackSplittingMonitor::doProfileX | ( | TH2 * | th2, |
MonitorElement * | me | ||
) | [private] |
void TrackSplittingMonitor::doProfileX | ( | MonitorElement * | th2m, |
MonitorElement * | me | ||
) | [private] |
void TrackSplittingMonitor::endJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 357 of file TrackSplittingMonitor.cc.
References conf_, dqmStore_, edm::ParameterSet::getParameter(), dumpDBToFile_GT_ttrig_cfg::outputFileName, DQMStore::save(), and DQMStore::showDirStructure().
{ bool outputMEsInRootFile = conf_.getParameter<bool>("OutputMEsInRootFile"); std::string outputFileName = conf_.getParameter<std::string>("OutputFileName"); if(outputMEsInRootFile){ dqmStore_->showDirStructure(); dqmStore_->save(outputFileName); } }
Definition at line 57 of file TrackSplittingMonitor.h.
Referenced by beginJob(), endJob(), and TrackSplittingMonitor().
Definition at line 62 of file TrackSplittingMonitor.h.
Referenced by analyze().
double TrackSplittingMonitor::d0Cut_ [private] |
Definition at line 72 of file TrackSplittingMonitor.h.
Referenced by analyze(), and beginJob().
Definition at line 98 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 84 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 105 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 91 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 93 of file TrackSplittingMonitor.h.
Referenced by analyze(), and beginJob().
Definition at line 79 of file TrackSplittingMonitor.h.
Referenced by analyze(), and beginJob().
Definition at line 100 of file TrackSplittingMonitor.h.
Referenced by analyze(), and beginJob().
Definition at line 86 of file TrackSplittingMonitor.h.
Referenced by analyze(), and beginJob().
Definition at line 94 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 80 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 101 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 87 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 95 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 81 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 102 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 88 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 97 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 83 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 104 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 90 of file TrackSplittingMonitor.h.
Referenced by beginJob().
DQMStore* TrackSplittingMonitor::dqmStore_ [private] |
Definition at line 56 of file TrackSplittingMonitor.h.
Referenced by beginJob(), endJob(), and TrackSplittingMonitor().
Definition at line 61 of file TrackSplittingMonitor.h.
Referenced by analyze().
Definition at line 96 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 82 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 103 of file TrackSplittingMonitor.h.
Referenced by beginJob().
Definition at line 89 of file TrackSplittingMonitor.h.
Referenced by beginJob().
double TrackSplittingMonitor::dzCut_ [private] |
Definition at line 73 of file TrackSplittingMonitor.h.
Referenced by analyze(), and beginJob().
std::string TrackSplittingMonitor::histname [private] |
Definition at line 54 of file TrackSplittingMonitor.h.
double TrackSplittingMonitor::norchiCut_ [private] |
Definition at line 75 of file TrackSplittingMonitor.h.
Referenced by analyze(), and beginJob().
int TrackSplittingMonitor::pixelHitsPerLeg_ [private] |
Definition at line 70 of file TrackSplittingMonitor.h.
Referenced by analyze(), and beginJob().
bool TrackSplittingMonitor::plotMuons_ [private] |
Definition at line 69 of file TrackSplittingMonitor.h.
Referenced by analyze(), and beginJob().
double TrackSplittingMonitor::ptCut_ [private] |
Definition at line 74 of file TrackSplittingMonitor.h.
Referenced by analyze(), and beginJob().
Definition at line 63 of file TrackSplittingMonitor.h.
Referenced by analyze().
Definition at line 66 of file TrackSplittingMonitor.h.
Referenced by analyze(), and beginJob().
Definition at line 65 of file TrackSplittingMonitor.h.
Referenced by analyze(), and beginJob().
Definition at line 59 of file TrackSplittingMonitor.h.
Referenced by analyze().
Definition at line 60 of file TrackSplittingMonitor.h.
Referenced by analyze().
int TrackSplittingMonitor::totalHitsPerLeg_ [private] |
Definition at line 71 of file TrackSplittingMonitor.h.
Referenced by analyze(), and beginJob().