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

chi2_vs_drj and pull_vs_pt validation plots #31400

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,7 +1,7 @@
import FWCore.ParameterSet.Config as cms

from SimTracker.TrackAssociatorProducers.quickTrackAssociatorByHits_cfi import *

from SimTracker.TrackAssociatorProducers.trackAssociatorByChi2_cfi import *
from SimTracker.TrackAssociation.trackingParticleRecoTrackAsssociation_cfi import *
import SimTracker.TrackAssociation.trackingParticleRecoTrackAsssociation_cfi
assoc2secStepTk = SimTracker.TrackAssociation.trackingParticleRecoTrackAsssociation_cfi.trackingParticleRecoTrackAsssociation.clone()
Expand Down
Expand Up @@ -81,7 +81,7 @@ SimToRecoCollection TrackAssociatorByChi2Impl::associateSimToReco(
const edm::RefToBaseVector<reco::Track>& tC, const edm::RefVector<TrackingParticleCollection>& tPCH) const {
const reco::BeamSpot& bs = *theBeamSpot;

SimToRecoCollection outputCollection;
SimToRecoCollection outputCollection(productGetter_);

int tpindex = 0;
for (auto tp = tPCH.begin(); tp != tPCH.end(); tp++, ++tpindex) {
Expand Down
Expand Up @@ -128,15 +128,17 @@ struct MTVHistoProducerAlgoForTrackerHistograms {
std::vector<METype> h_assochi2, h_assochi2_prob;

//chi2 and # lost hits vs eta: to be used with doProfileX
std::vector<METype> chi2_vs_eta, chi2_vs_pt, nlosthits_vs_eta;
std::vector<METype> assoc_chi2_vs_eta, assoc_chi2_vs_pt, assoc_chi2prob_vs_eta, assoc_chi2prob_vs_pt;
std::vector<METype> chi2_vs_eta, chi2_vs_pt, chi2_vs_drj, nlosthits_vs_eta;
std::vector<METype> assoc_chi2_vs_eta, assoc_chi2_vs_pt, assoc_chi2_vs_drj, assoc_chi2prob_vs_eta,
assoc_chi2prob_vs_pt, assoc_chi2prob_vs_drj;

//resolution of track params: to be used with fitslicesytool
std::vector<METype> dxyres_vs_eta, ptres_vs_eta, dzres_vs_eta, phires_vs_eta, cotThetares_vs_eta;
std::vector<METype> dxyres_vs_pt, ptres_vs_pt, dzres_vs_pt, phires_vs_pt, cotThetares_vs_pt;

//pulls of track params vs eta: to be used with fitslicesytool
std::vector<METype> dxypull_vs_eta, ptpull_vs_eta, dzpull_vs_eta, phipull_vs_eta, thetapull_vs_eta;
std::vector<METype> dxypull_vs_pt, ptpull_vs_pt, dzpull_vs_pt, phipull_vs_pt, thetapull_vs_pt;
std::vector<METype> ptpull_vs_phi, phipull_vs_phi, thetapull_vs_phi;
};

Expand Down
Expand Up @@ -82,8 +82,8 @@
#
# dR_jet
mindrj = cms.double(0.001),
maxdrj = cms.double(0.1),
nintdrj = cms.int32(50),
maxdrj = cms.double(0.5),
nintdrj = cms.int32(250),
#
# chi2/ndof
minChi2 = cms.double(0),
Expand Down
5 changes: 5 additions & 0 deletions Validation/RecoTrack/python/PostProcessorTracker_cfi.py
Expand Up @@ -173,25 +173,30 @@ def _addNoFlow(module):
"cotThetares_vs_eta '#sigma(cot(#theta)) vs #eta' cotThetares_vs_eta",
"cotThetares_vs_pt '#sigma(cot(#theta)) vs p_{T}' cotThetares_vs_pt",
"h_dxypulleta 'd_{xy} Pull vs #eta' dxypull_vs_eta",
"h_dxypullpt 'd_{xy} Pull vs p_{T}' dxypull_vs_pt",
"dxyres_vs_eta '#sigma(d_{xy}) vs #eta' dxyres_vs_eta",
"dxyres_vs_phi '#sigma(d_{xy}) vs #phi' dxyres_vs_phi",
"dxyres_vs_pt '#sigma(d_{xy}) vs p_{T}' dxyres_vs_pt",
"h_dzpulleta 'd_{z} Pull vs #eta' dzpull_vs_eta",
"h_dzpullpt 'd_{z} Pull vs p_{T}' dzpull_vs_pt",
"dzres_vs_eta '#sigma(d_{z}) vs #eta' dzres_vs_eta",
"dzres_vs_phi '#sigma(d_{z}) vs #phi' dzres_vs_phi",
"dzres_vs_pt '#sigma(d_{z}) vs p_{T}' dzres_vs_pt",
"etares_vs_eta '#sigma(#eta) vs #eta' etares_vs_eta",
"h_phipulleta '#phi Pull vs #eta' phipull_vs_eta",
"h_phipullpt '#phi Pull vs p_{T}' phipull_vs_pt",
"h_phipullphi '#phi Pull vs #phi' phipull_vs_phi",
"phires_vs_eta '#sigma(#phi) vs #eta' phires_vs_eta",
"phires_vs_phi '#sigma(#phi) vs #phi' phires_vs_phi",
"phires_vs_pt '#sigma(#phi) vs p_{T}' phires_vs_pt",
"h_ptpulleta 'p_{T} Pull vs #eta' ptpull_vs_eta",
"h_ptpullpt 'p_{T} Pull vs p_{T}' ptpull_vs_pt",
"h_ptpullphi 'p_{T} Pull vs #phi' ptpull_vs_phi",
"ptres_vs_eta '#sigma(p_{T}) vs #eta' ptres_vs_eta",
"ptres_vs_phi '#sigma(p_{T}) vs #phi' ptres_vs_phi",
"ptres_vs_pt '#sigma(p_{T}) vs p_{T}' ptres_vs_pt",
"h_thetapulleta '#theta Pull vs #eta' thetapull_vs_eta",
"h_thetapullpt '#theta Pull vs p_{T}' thetapull_vs_pt",
"h_thetapullphi '#theta Pull vs #phi' thetapull_vs_phi"
),
cumulativeDists = cms.untracked.vstring(
Expand Down
15 changes: 14 additions & 1 deletion Validation/RecoTrack/python/TrackValidation_cff.py
@@ -1,7 +1,7 @@
from __future__ import absolute_import
import FWCore.ParameterSet.Config as cms

import SimTracker.TrackAssociatorProducers.trackAssociatorByChi2_cfi
from SimTracker.TrackAssociatorProducers.trackAssociatorByChi2_cfi import *
from SimTracker.TrackAssociatorProducers.quickTrackAssociatorByHits_cfi import *
from SimTracker.TrackAssociation.trackingParticleRecoTrackAsssociation_cfi import *
import Validation.RecoTrack.MultiTrackValidator_cfi
Expand Down Expand Up @@ -353,6 +353,11 @@ def _getMVASelectors(postfix):
ptMin = 0,
)

#ByChi2 association (for jetCore usage, not used by default)
MTVTrackAssociationByChi2 = trackingParticleRecoTrackAsssociation.clone(
associator = cms.InputTag('trackAssociatorByChi2')
)

# Select jets for JetCore tracking
highPtJets = cms.EDFilter("CandPtrSelector", src = cms.InputTag("ak4CaloJets"), cut = cms.string("pt()>1000"))
highPtJetsForTrk = highPtJetsForTrk = highPtJets.clone(src = "ak4CaloJetsForTrk")
Expand All @@ -364,6 +369,7 @@ def _getMVASelectors(postfix):
trackValidator = Validation.RecoTrack.MultiTrackValidator_cfi.multiTrackValidator.clone(
useLogPt = cms.untracked.bool(True),
dodEdxPlots = True,
# associators=cms.untracked.VInputTag('MTVTrackAssociationByChi2'), #uncomment for byChi2 assoc. for jetcore studies (1/5)
doPVAssociationPlots = True
#,minpT = cms.double(-1)
#,maxpT = cms.double(3)
Expand All @@ -387,6 +393,7 @@ def _getMVASelectors(postfix):
locals()["_generalTracksHp"+_postfix],
"generalTracksPt09",
"cutsRecoTracksBtvLike",
"cutsRecoTracksJetCoreRegionalStepByOriginalAlgo",
]
)
_setForEra(trackValidator.histoProducerAlgoBlock, _eraName, _era, seedingLayerSets=locals()["_seedingLayerSets"+_postfix])
Expand Down Expand Up @@ -532,6 +539,9 @@ def _getMVASelectors(postfix):
trackValidatorBuilding = _trackValidatorSeedingBuilding.clone(
dirName = "Tracking/TrackBuilding/",
doMVAPlots = True,
doResolutionPlotsForLabels = ['jetCoreRegionalStepTracks'],
# associators = ["trackAssociatorByChi2"], #uncomment for byChi2 assoc. for jetcore studies (2/5)
# UseAssociators = True, #uncomment for byChi2 assoc. for jetcore studies (3/5)
)
trackValidatorBuildingPreSplitting = trackValidatorBuilding.clone(
associators = ["quickTrackAssociatorByHitsPreSplitting"],
Expand Down Expand Up @@ -648,6 +658,8 @@ def _uniqueFirstLayers(layerList):
tracksValidationTruth = cms.Task(
tpClusterProducer,
tpClusterProducerPreSplitting,
# trackAssociatorByChi2, #uncomment for byChi2 assoc. for jetcore studies (4/5)
# MTVTrackAssociationByChi2, #uncomment for byChi2 assoc. for jetcore studies (5/5)
quickTrackAssociatorByHits,
quickTrackAssociatorByHitsPreSplitting,
trackingParticleRecoTrackAsssociation,
Expand Down Expand Up @@ -829,6 +841,7 @@ def _uniqueFirstLayers(layerList):
dirName = "Tracking/TrackSeeding/",
label = _seedSelectors,
doSeedPlots = True,
doResolutionPlotsForLabels = [ "seedTracksjetCoreRegionalStepSeeds",]
)
trackValidatorSeedingPreSplittingTrackingOnly = trackValidatorSeedingTrackingOnly.clone(
associators = ["quickTrackAssociatorByHitsPreSplitting"],
Expand Down
1 change: 1 addition & 0 deletions Validation/RecoTrack/python/associators_cff.py
Expand Up @@ -2,6 +2,7 @@

#### TrackAssociation
import SimTracker.TrackAssociatorProducers.quickTrackAssociatorByHits_cfi
import SimTracker.TrackAssociatorProducers.trackAssociatorByChi2_cfi
import SimTracker.TrackAssociatorProducers.trackAssociatorByPosition_cfi
from SimTracker.TrackerHitAssociation.tpClusterProducer_cfi import tpClusterProducer as _tpClusterProducer
from SimTracker.TrackAssociation.trackingParticleRecoTrackAsssociation_cfi import trackingParticleRecoTrackAsssociation as _trackingParticleRecoTrackAsssociation
Expand Down
10 changes: 6 additions & 4 deletions Validation/RecoTrack/python/plotting/trackingPlots.py
Expand Up @@ -41,6 +41,7 @@
_maxPU = [20, 50, 65, 80, 100, 150, 200, 250]
_minMaxTracks = [0, 200, 500, 1000, 1500, 2000]
_minMaxMVA = [-1.025, -0.5, 0, 0.5, 1.025]
_maxDRJ = 0.1

def _minMaxResidual(ma):
return ([-x for x in ma], ma)
Expand Down Expand Up @@ -216,7 +217,7 @@ def _makeMVAPlots(num, hp=False):
)
_effandfakeDeltaRPU = PlotGroup("effandfakeDeltaRPU",
_makeEffFakeDupPlots("dr" , "#DeltaR", effopts=dict(xtitle="TP min #DeltaR"), fakeopts=dict(xtitle="track min #DeltaR"), common=dict(xlog=True)) +
_makeEffFakeDupPlots("drj" , "#DeltaR(track, jet)", effopts=dict(xtitle="#DeltaR(TP, jet)", ytitle="efficiency vs #DeltaR(TP, jet"), fakeopts=dict(xtitle="#DeltaR(track, jet)"), common=dict(xlog=True))+
_makeEffFakeDupPlots("drj" , "#DeltaR(track, jet)", effopts=dict(xtitle="#DeltaR(TP, jet)", ytitle="efficiency vs #DeltaR(TP, jet"), fakeopts=dict(xtitle="#DeltaR(track, jet)"), common=dict(xlog=True, xmax=_maxDRJ))+
_makeEffFakeDupPlots("pu" , "PU" , common=dict(xtitle="Pileup", xmin=_minPU, xmax=_maxPU)),
legendDy=_legendDy_4rows
)
Expand Down Expand Up @@ -262,7 +263,7 @@ def _makeMVAPlots(num, hp=False):
)
_dupandfakeDeltaRPU = PlotGroup("dupandfakeDeltaRPU",
_makeFakeDupPileupPlots("dr" , "#DeltaR", xquantity="min #DeltaR", common=dict(xlog=True)) +
_makeFakeDupPileupPlots("drj" , "#DeltaR(track, jet)", xtitle="#DeltaR(track, jet)", common=dict(xlog=True)) +
_makeFakeDupPileupPlots("drj" , "#DeltaR(track, jet)", xtitle="#DeltaR(track, jet)", common=dict(xlog=True, xmax=_maxDRJ)) +
_makeFakeDupPileupPlots("pu" , "PU" , xtitle="Pileup", common=dict(xmin=_minPU, xmax=_maxPU)),
ncols=3
)
Expand Down Expand Up @@ -373,6 +374,7 @@ def _makeMVAPlots(num, hp=False):
fallback={"name": "chi2_vs_eta", "profileX": True}),
Plot("ptres_vs_eta_Mean", scale=100, title="", xtitle="TP #eta (PCA to beamline)", ytitle="< #delta p_{T} / p_{T} > (%)", ymin=_minResidualPt, ymax=_maxResidualPt),
Plot("chi2mean_vs_pt", title="", xtitle="p_{T}", ytitle="< #chi^{2} / ndf >", ymin=[0, 0.5], ymax=[2, 2.5, 3, 5], xlog=True, fallback={"name": "chi2_vs_pt", "profileX": True}),
Plot("chi2mean_vs_drj", title="", xtitle="#DeltaR(track, jet)", ytitle="< #chi^{2} / ndf >", ymin=[0, 0.5], ymax=[2, 2.5, 3, 5], xlog=True, xmax=_maxDRJ, fallback={"name": "chi2_vs_drj", "profileX": True}),
Plot("ptres_vs_pt_Mean", title="", xtitle="p_{T}", ytitle="< #delta p_{T}/p_{T} > (%)", scale=100, ymin=_minResidualPt, ymax=_maxResidualPt,xlog=True)
])
_common = {"stat": True, "fit": True, "normalizeToUnitArea": True, "drawStyle": "hist", "drawCommand": "", "xmin": -10, "xmax": 10, "ylog": True, "ymin": 5e-5, "ymax": [0.01, 0.05, 0.1, 0.2, 0.5, 0.8, 1.025], "ratioUncertainty": False}
Expand Down Expand Up @@ -441,7 +443,7 @@ def _makeMVAPlots(num, hp=False):
)
_extDistDeltaR = PlotGroup("distDeltaR",
_makeDistPlots("dr" , "min #DeltaR", common=dict(xlog=True)) +
_makeDistPlots("drj" , "#DeltaR(track, jet)", common=dict(xlog=True)),
_makeDistPlots("drj" , "#DeltaR(track, jet)", common=dict(xlog=True, xmax=_maxDRJ)),
ncols=2, legendDy=_legendDy_2rows,
)
_extDistSeedingPlots = _makeDistPlots("seedingLayerSet", "seeding layers", common=dict(xtitle="", **_seedingLayerSet_common))
Expand Down Expand Up @@ -512,7 +514,7 @@ def _makeMVAPlots(num, hp=False):
)
_extDistSimDeltaR = PlotGroup("distsimDeltaR",
_makeDistSimPlots("dr" , "min #DeltaR", common=dict(xlog=True)) +
_makeDistSimPlots("drj" , "#DeltaR(TP, jet)", common=dict(xlog=True)),
_makeDistSimPlots("drj" , "#DeltaR(TP, jet)", common=dict(xlog=True, xmax=_maxDRJ)),
ncols=2, legendDy=_legendDy_2rows,
)

Expand Down
45 changes: 44 additions & 1 deletion Validation/RecoTrack/src/MTVHistoProducerAlgoForTracker.cc
Expand Up @@ -1099,6 +1099,8 @@ void MTVHistoProducerAlgoForTracker::bookRecoHistos(DQMStore::IBooker& ibook,
ibook.bookProfile("chi2mean_vs_phi", "mean #chi^{2} vs #phi", nintPhi, minPhi, maxPhi, 200, 0, 20, " "));
histograms.chi2_vs_pt.push_back(
makeProfileIfLogX(ibook, useLogPt, "chi2mean_vs_pt", "mean #chi^{2} vs p_{T}", nintPt, minPt, maxPt, 0, 20));
histograms.chi2_vs_drj.push_back(makeProfileIfLogX(
ibook, true, "chi2mean_vs_drj", "mean #chi^{2} vs dR(track,jet)", nintdrj, log10(mindrj), log10(maxdrj), 0, 20));

histograms.assoc_chi2_vs_eta.push_back(
ibook.bookProfile("assoc_chi2mean", "mean #chi^{2} vs #eta", nintEta, minEta, maxEta, 200, 0., 20., " "));
Expand All @@ -1107,7 +1109,25 @@ void MTVHistoProducerAlgoForTracker::bookRecoHistos(DQMStore::IBooker& ibook,
histograms.assoc_chi2_vs_pt.push_back(makeProfileIfLogX(
ibook, useLogPt, "assoc_chi2mean_vs_pt", "mean #chi^{2} vs p_{T}", nintPt, minPt, maxPt, 0., 20.));
histograms.assoc_chi2prob_vs_pt.push_back(makeProfileIfLogX(
ibook, useLogPt, "assoc_chi2prob_vs_pt", "mean #chi^{2} probability vs p_{T}", nintPt, minPt, maxPt, 0., 20.));
Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

here the only change in the plots I've done (max Y limit from 20 to 1, because it's a probability. I thought it was a typo from another PR copying from the line above). If I am wrong I change it back immediately.

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I agree with this update, thanks

ibook, useLogPt, "assoc_chi2prob_vs_pt", "mean #chi^{2} probability vs p_{T}", nintPt, minPt, maxPt, 0., 1.));
histograms.assoc_chi2_vs_drj.push_back(makeProfileIfLogX(ibook,
true,
"assoc_chi2mean_vs_drj",
"mean #chi^{2} vs dR(track,jet)",
nintdrj,
log10(mindrj),
log10(maxdrj),
0.,
20));
histograms.assoc_chi2prob_vs_drj.push_back(makeProfileIfLogX(ibook,
true,
"assoc_chi2prob_vs_drj",
"mean #chi^{2} probability vs dR(track,jet)",
nintdrj,
log10(mindrj),
log10(maxdrj),
0.,
1.));

histograms.nhits_vs_eta.push_back(
ibook.bookProfile("hits_eta", "mean hits vs eta", nintEta, minEta, maxEta, nintHit, minHit, maxHit, " "));
Expand Down Expand Up @@ -1387,6 +1407,16 @@ void MTVHistoProducerAlgoForTracker::bookRecoHistos(DQMStore::IBooker& ibook,
histograms.phipull_vs_eta, false, "phipull_vs_eta", "phipull_vs_eta", nintEta, minEta, maxEta, 100, -10, 10);
bookResolutionPlots2D(
histograms.thetapull_vs_eta, false, "thetapull_vs_eta", "thetapull_vs_eta", nintEta, minEta, maxEta, 100, -10, 10);
bookResolutionPlots2D(
histograms.dxypull_vs_pt, useLogPt, "dxypull_vs_pt", "dxypull_vs_pt", nintPt, minPt, maxPt, 100, -10, 10);
bookResolutionPlots2D(
histograms.ptpull_vs_pt, useLogPt, "ptpull_vs_pt", "ptpull_vs_pt", nintPt, minPt, maxPt, 100, -10, 10);
bookResolutionPlots2D(
histograms.dzpull_vs_pt, useLogPt, "dzpull_vs_pt", "dzpull_vs_pt", nintPt, minPt, maxPt, 100, -10, 10);
bookResolutionPlots2D(
histograms.phipull_vs_pt, useLogPt, "phipull_vs_pt", "phipull_vs_pt", nintPt, minPt, maxPt, 100, -10, 10);
bookResolutionPlots2D(
histograms.thetapull_vs_pt, useLogPt, "thetapull_vs_pt", "thetapull_vs_pt", nintPt, minPt, maxPt, 100, -10, 10);

// histograms.h_ptshiftetamean.push_back( ibook.book1D("h_ptshifteta_Mean","<#deltapT/pT>[%] vs #eta",nintEta,minEta,maxEta) );

Expand Down Expand Up @@ -2026,6 +2056,8 @@ void MTVHistoProducerAlgoForTracker::fill_generic_recoTrack_histos(const Histogr
histograms.h_recozpos[count]->Fill(vertz);
histograms.h_recodr[count]->Fill(dR);
histograms.h_recodrj[count]->Fill(dRJet);
if (dRJet <= 99999) //dRJet can be set to numeric_limits max^2, this is a protection
histograms.chi2_vs_drj[count]->Fill(dRJet, chi2);
if (fillSeedingLayerSets)
histograms.h_reco_seedingLayerSet[count]->Fill(seedingLayerSetBin);
if (pvPosition) {
Expand Down Expand Up @@ -2089,6 +2121,10 @@ void MTVHistoProducerAlgoForTracker::fill_generic_recoTrack_histos(const Histogr
histograms.assoc_chi2prob_vs_eta[count]->Fill(eta, chi2prob);
histograms.assoc_chi2_vs_pt[count]->Fill(pt, chi2);
histograms.assoc_chi2prob_vs_pt[count]->Fill(pt, chi2prob);
if (dRJet <= 99999) { //dRJet can be set to numeric_limits max^2, this is a protection
histograms.assoc_chi2_vs_drj[count]->Fill(dRJet, chi2);
histograms.assoc_chi2prob_vs_drj[count]->Fill(dRJet, chi2prob);
}
histograms.h_assoc2vertpos[count]->Fill(vertxy);
histograms.h_assoc2zpos[count]->Fill(vertz);
histograms.h_assoc2dr[count]->Fill(dR);
Expand Down Expand Up @@ -2481,6 +2517,13 @@ void MTVHistoProducerAlgoForTracker::fill_ResoAndPull_recoTrack_histos(const His
histograms.phipull_vs_eta[count]->Fill(etaSim, phiPull);
histograms.thetapull_vs_eta[count]->Fill(etaSim, thetaPull);

//pulls of track params vs pt: fill 2D histos
histograms.dxypull_vs_pt[count]->Fill(ptSim, dxyPull);
histograms.ptpull_vs_pt[count]->Fill(ptSim, ptres / ptError);
histograms.dzpull_vs_pt[count]->Fill(ptSim, dzPull);
histograms.phipull_vs_pt[count]->Fill(ptSim, phiPull);
histograms.thetapull_vs_pt[count]->Fill(ptSim, thetaPull);

//plots vs phi
histograms.nhits_vs_phi[count]->Fill(phiRec, track.numberOfValidHits());
histograms.chi2_vs_phi[count]->Fill(phiRec, track.normalizedChi2());
Expand Down
1 change: 1 addition & 0 deletions Validation/RecoTrack/test/MultiTrackValidator_cfg.py
Expand Up @@ -96,6 +96,7 @@
### validation-specific includes
#process.load("SimTracker.TrackAssociatorProducers.trackAssociatorByHits_cfi")
process.load("SimTracker.TrackAssociatorProducers.quickTrackAssociatorByHits_cfi")
process.load("SimTracker.TrackAssociatorProducers.trackAssociatorByChi2_cfi")
process.load("SimTracker.TrackAssociation.trackingParticleRecoTrackAsssociation_cfi")
process.load("Validation.RecoTrack.cuts_cff")
process.load("Validation.RecoTrack.MultiTrackValidator_cff")
Expand Down