Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add impact parameter significance histograms to tracking DQM #8410

Merged
merged 1 commit into from Mar 20, 2015
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
1 change: 1 addition & 0 deletions DQM/TrackingMonitor/BuildFile.xml
Expand Up @@ -10,6 +10,7 @@
<use name="TrackingTools/TransientTrackingRecHit"/>
<use name="TrackingTools/TransientTrack"/>
<use name="TrackingTools/DetLayers"/>
<use name="TrackingTools/IPTools"/>
<use name="TrackPropagation/SteppingHelixPropagator"/>
<use name="DataFormats/MuonReco"/>
<use name="Geometry/RPCGeometry"/>
Expand Down
13 changes: 12 additions & 1 deletion DQM/TrackingMonitor/interface/TrackAnalyzer.h
Expand Up @@ -94,6 +94,10 @@ class TrackAnalyzer

//For HI Plots
bool doHIPlots_;

// IP significance plots
bool doSIPPlots_;

std::string qualityString_;

struct TkParameterMEs {
Expand Down Expand Up @@ -267,7 +271,14 @@ class TrackAnalyzer
MonitorElement* dNdPt_HighPurity;
MonitorElement* NhitVsEta_HighPurity;
MonitorElement* NhitVsPhi_HighPurity;


// IP significance plots
MonitorElement *sipDxyToBS;
MonitorElement *sipDzToBS;
MonitorElement *sip3dToPV;
MonitorElement *sip2dToPV;
MonitorElement *sipDxyToPV;
MonitorElement *sipDzToPV;

struct TkRecHitsPerSubDetMEs {
MonitorElement* NumberOfRecHitsPerTrack;
Expand Down
1 change: 1 addition & 0 deletions DQM/TrackingMonitor/python/TrackingMonitor_cfi.py
Expand Up @@ -57,6 +57,7 @@
doDCAPlots = cms.bool(False),
doDCAwrtPVPlots = cms.bool(False),
doDCAwrt000Plots = cms.bool(False),
doSIPPlots = cms.bool(False),
doGeneralPropertiesPlots = cms.bool(False),
doHitPropertiesPlots = cms.bool(False),
# doGoodTrackPlots = cms.bool(False),
Expand Down
77 changes: 73 additions & 4 deletions DQM/TrackingMonitor/src/TrackAnalyzer.cc
Expand Up @@ -9,6 +9,7 @@
#include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
#include "TrackingTools/Records/interface/TransientTrackRecord.h"
#include "TrackingTools/TransientTrack/interface/TransientTrack.h"
#include "TrackingTools/IPTools/interface/IPTools.h"
#include "MagneticField/Engine/interface/MagneticField.h"

#include "FWCore/MessageLogger/interface/MessageLogger.h"
Expand Down Expand Up @@ -38,6 +39,7 @@ TrackAnalyzer::TrackAnalyzer(const edm::ParameterSet& iConfig)
, doLumiAnalysis_ ( conf_.getParameter<bool>("doLumiAnalysis") )
, doTestPlots_ ( conf_.getParameter<bool>("doTestPlots") )
, doHIPlots_ ( conf_.getParameter<bool>("doHIPlots") )
, doSIPPlots_ ( conf_.getParameter<bool>("doSIPPlots") )
, qualityString_ ( conf_.getParameter<std::string>("qualityString"))
{
initHistos();
Expand Down Expand Up @@ -118,6 +120,13 @@ void TrackAnalyzer::initHistos()
NhitVsEta_HighPurity = NULL;
NhitVsPhi_HighPurity = NULL;

// IP significance
sipDxyToBS = NULL;
sipDzToBS = NULL;
sip3dToPV = NULL;
sip2dToPV = NULL;
sipDxyToPV = NULL;
sipDzToPV = NULL;

}

Expand Down Expand Up @@ -622,6 +631,46 @@ void TrackAnalyzer::bookHistosForBeamSpot(DQMStore::IBooker & ibooker) {
}
}


if (doSIPPlots_ || doAllPlots_) {
const double sipBins = 200;
const double sipMin = -20;
const double sipMax = 20;

ibooker.setCurrentFolder(TopFolder_+"/GeneralProperties");

// SIP wrt. beamspot
histname = "SIPDxyToBS_";
sipDxyToBS = ibooker.book1D(histname+CategoryName, histname+CategoryName, sipBins, sipMin, sipMax);
sipDxyToBS->setAxisTitle("Track dxy significance wrt beam spot",1);
sipDxyToBS->setAxisTitle("Number of Tracks",2);

histname = "SIPDzToBS_";
sipDzToBS = ibooker.book1D(histname+CategoryName, histname+CategoryName, sipBins, sipMin, sipMax);
sipDzToBS->setAxisTitle("Track dz significance wrt beam spot",1);
sipDzToBS->setAxisTitle("Number of Tracks",2);

// SIP wrt. vertex
histname = "SIP3DToPV_";
sip3dToPV = ibooker.book1D(histname+CategoryName, histname+CategoryName, sipBins, sipMin, sipMax);
sip3dToPV->setAxisTitle("3D IP significance wrt primary vertex",1);
sip3dToPV->setAxisTitle("Number of Tracks",2);

histname = "SIP2DToPV_";
sip2dToPV = ibooker.book1D(histname+CategoryName, histname+CategoryName, sipBins, sipMin, sipMax);
sip2dToPV->setAxisTitle("2D IP significance wrt primary vertex",1);
sip2dToPV->setAxisTitle("Number of Tracks",2);

histname = "SIPDxyToPV_";
sipDxyToPV = ibooker.book1D(histname+CategoryName, histname+CategoryName, sipBins, sipMin, sipMax);
sipDxyToPV->setAxisTitle("Track dxy significance wrt primary vertex",1);
sipDxyToPV->setAxisTitle("Number of Tracks",2);

histname = "SIPDzToPV_";
sipDzToPV = ibooker.book1D(histname+CategoryName, histname+CategoryName, sipBins, sipMin, sipMax);
sipDzToPV->setAxisTitle("Track dz significance wrt primary vertex",1);
sipDzToPV->setAxisTitle("Number of Tracks",2);
}
}

// -- Analyse
Expand Down Expand Up @@ -687,11 +736,11 @@ void TrackAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSe
Chi2oNDF_lumiFlag -> Fill(chi2oNDF);
}

if(doDCAPlots_ || doBSPlots_ || doAllPlots_) {
if(doDCAPlots_ || doBSPlots_ || doSIPPlots_ || doAllPlots_) {

edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
iEvent.getByToken(beamSpotToken_,recoBeamSpotHandle);
reco::BeamSpot bs = *recoBeamSpotHandle;
const reco::BeamSpot& bs = *recoBeamSpotHandle;

DistanceOfClosestApproachToBS -> Fill(track.dxy(bs.position()));
DistanceOfClosestApproachToBSVsPhi -> Fill(track.phi(), track.dxy(bs.position()));
Expand All @@ -704,13 +753,18 @@ void TrackAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSe
TESTDistanceOfClosestApproachToBS -> Fill(track.dxy(bs.position(track.vz())));
TESTDistanceOfClosestApproachToBSVsPhi -> Fill(track.phi(), track.dxy(bs.position(track.vz())));
}

if(doSIPPlots_) {
sipDxyToBS->Fill(track.dxy(bs.position())/track.dxyError());
sipDzToBS->Fill(track.dz(bs.position())/track.dzError());
}
}

if(doDCAPlots_ || doPVPlots_ || doAllPlots_) {
if(doDCAPlots_ || doPVPlots_ || doSIPPlots_ || doAllPlots_) {
edm::Handle<reco::VertexCollection> recoPrimaryVerticesHandle;
iEvent.getByToken(pvToken_,recoPrimaryVerticesHandle);
if (recoPrimaryVerticesHandle->size() > 0) {
reco::Vertex pv = recoPrimaryVerticesHandle->at(0);
const reco::Vertex& pv = (*recoPrimaryVerticesHandle)[0];


//////////////////
Expand Down Expand Up @@ -742,6 +796,21 @@ void TrackAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSe
DistanceOfClosestApproachToPVVsPhi -> Fill(track.phi(), track.dxy(pv.position()));
xPointOfClosestApproachVsZ0wrtPV -> Fill(track.dz(pv.position()),(track.vx()-pv.position().x()));
yPointOfClosestApproachVsZ0wrtPV -> Fill(track.dz(pv.position()),(track.vy()-pv.position().y()));


if(doSIPPlots_) {
edm::ESHandle<TransientTrackBuilder> theB;
iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
reco::TransientTrack transTrack = theB->build(track);

GlobalVector dir(track.px(), track.py(), track.pz());
std::pair<bool, Measurement1D> ip3d = IPTools::signedImpactParameter3D(transTrack, dir, pv);
std::pair<bool, Measurement1D> ip2d = IPTools::signedTransverseImpactParameter(transTrack, dir, pv);
sip3dToPV->Fill(ip3d.second.value() / ip3d.second.error());
sip2dToPV->Fill(ip2d.second.value() / ip2d.second.error());
sipDxyToPV->Fill(track.dxy(pv.position())/track.dxyError());
sipDzToPV->Fill(track.dz(pv.position())/track.dzError());
}
}
}

Expand Down
Expand Up @@ -23,6 +23,7 @@
locals()[label].doDCAPlots = doPlotsPCA[tracks]
locals()[label].doDCAwrtPVPlots = doPlotsPCA[tracks]
locals()[label].doDCAwrt000Plots = doPlotsPCA[tracks]
locals()[label].doSIPPlots = doPlotsPCA[tracks]
locals()[label].numCut = numCutString[tracks]
locals()[label].denCut = denCutString[tracks]
locals()[label].setLabel(label)
Expand All @@ -40,6 +41,7 @@
locals()[label].doDCAPlots = doPlotsPCA[tracks]
locals()[label].doDCAwrtPVPlots = doPlotsPCA[tracks]
locals()[label].doDCAwrt000Plots = doPlotsPCA[tracks]
locals()[label].doSIPPlots = doPlotsPCA[tracks]
locals()[label].numCut = numCutString[tracks]
locals()[label].denCut = denCutString[tracks]
locals()[label].setLabel(label)
Expand Down