From f61b144d5bc13c963299934b0eefd55972f2e853 Mon Sep 17 00:00:00 2001 From: mmusich Date: Tue, 2 Jan 2024 16:46:46 +0100 Subject: [PATCH] miscellaneous improvements to SingleLongTrackProducer. Also use earlyMuons collections in order to comply with the trackingOnly DQM sequence --- .../python/shortTrackResolution_cff.py | 27 ++++++---- .../plugins/SingleLongTrackProducer.cc | 51 ++++++++++++++----- .../python/SingleLongTrackProducer_cfi.py | 2 +- 3 files changed, 56 insertions(+), 24 deletions(-) diff --git a/DQM/TrackingMonitorSource/python/shortTrackResolution_cff.py b/DQM/TrackingMonitorSource/python/shortTrackResolution_cff.py index e38489adb17f9..b07c24d88cb99 100644 --- a/DQM/TrackingMonitorSource/python/shortTrackResolution_cff.py +++ b/DQM/TrackingMonitorSource/python/shortTrackResolution_cff.py @@ -5,32 +5,32 @@ from RecoTracker.FinalTrackSelectors.trackerTrackHitFilter_cfi import trackerTrackHitFilter as _trackerTrackHitFilter ShortTrackCandidates = _trackerTrackHitFilter.clone(src = "SingleLongTrackProducer", - truncateTracks = True, - replaceWithInactiveHits = True, - rejectBadStoNHits = True, - usePixelQualityFlag = True) + truncateTracks = True, + replaceWithInactiveHits = True, + rejectBadStoNHits = True, + usePixelQualityFlag = True) from Configuration.Eras.Modifier_phase2_tracker_cff import phase2_tracker phase2_tracker.toModify(ShortTrackCandidates, isPhase2 = True) ShortTrackCandidates3 = ShortTrackCandidates.clone(minimumHits = 3, - layersRemaining = 3) + layersRemaining = 3) ShortTrackCandidates4 = ShortTrackCandidates.clone(minimumHits = 4, - layersRemaining = 4) + layersRemaining = 4) ShortTrackCandidates5 = ShortTrackCandidates.clone(minimumHits = 5, - layersRemaining = 5) + layersRemaining = 5) ShortTrackCandidates6 = ShortTrackCandidates.clone(minimumHits = 6, - layersRemaining = 6) + layersRemaining = 6) ShortTrackCandidates7 = ShortTrackCandidates.clone(minimumHits = 7, - layersRemaining = 7) + layersRemaining = 7) ShortTrackCandidates8 = ShortTrackCandidates.clone(minimumHits = 8, - layersRemaining = 8) + layersRemaining = 8) import RecoTracker.TrackProducer.CTFFinalFitWithMaterial_cff RefittedShortTracks = RecoTracker.TrackProducer.CTFFinalFitWithMaterial_cff.ctfWithMaterialTracks.clone(src = 'ShortTrackCandidates') @@ -51,7 +51,12 @@ maxTracksPtInput = 99999.9, maxDrInput = 0.01, tracksInputTag = "SingleLongTrackProducer", - tracksRerecoInputTag = ["RefittedShortTracks3","RefittedShortTracks4","RefittedShortTracks5","RefittedShortTracks6","RefittedShortTracks7","RefittedShortTracks8"]) + tracksRerecoInputTag = ["RefittedShortTracks3", + "RefittedShortTracks4", + "RefittedShortTracks5", + "RefittedShortTracks6", + "RefittedShortTracks7", + "RefittedShortTracks8"]) shortTrackResolution3to8 = cms.Sequence(SingleLongTrackProducer * ShortTrackCandidates3 * diff --git a/RecoTracker/FinalTrackSelectors/plugins/SingleLongTrackProducer.cc b/RecoTracker/FinalTrackSelectors/plugins/SingleLongTrackProducer.cc index 7af6c8c97145b..a9b8141f423ac 100644 --- a/RecoTracker/FinalTrackSelectors/plugins/SingleLongTrackProducer.cc +++ b/RecoTracker/FinalTrackSelectors/plugins/SingleLongTrackProducer.cc @@ -54,15 +54,35 @@ SingleLongTrackProducer::SingleLongTrackProducer(const edm::ParameterSet &iConfi void SingleLongTrackProducer::produce(edm::Event &iEvent, const edm::EventSetup &iSetup) { using namespace edm; + // register output collection: + std::unique_ptr goodTracks(new reco::TrackCollection); + // register input collections: const auto &tracks = iEvent.getHandle(tracksToken); + if (!tracks.isValid()) { + edm::LogError("SingleLongTrackProducer") + << "Input track collection is not valid.\n Returning empty output track collection."; + iEvent.put(std::move(goodTracks), ""); + return; + } + const auto &muons = iEvent.getHandle(muonsToken); + if (!muons.isValid()) { + edm::LogError("SingleLongTrackProducer") + << "Input muon collection is not valid.\n Returning empty output track collection."; + iEvent.put(std::move(goodTracks), ""); + return; + } - const auto &vtx_h = iEvent.getHandle(PrimVtxToken); - const reco::Vertex vtx = vtx_h->at(0); + const auto &vertices = iEvent.getHandle(PrimVtxToken); + if (!vertices.isValid()) { + edm::LogError("SingleLongTrackProducer") + << "Input vertex collection is not valid.\n Returning empty output track collection."; + iEvent.put(std::move(goodTracks), ""); + return; + } - // register output collection: - std::unique_ptr goodTracks(new reco::TrackCollection); + const reco::Vertex vtx = vertices->at(0); // Preselection of long quality tracks std::vector selTracks; @@ -80,7 +100,7 @@ void SingleLongTrackProducer::produce(edm::Event &iEvent, const edm::EventSetup if (track.pt() < minPt) continue; - if (abs(track.eta()) > maxEta) + if (std::abs(track.eta()) > maxEta) continue; if (hitpattern.trackerLayersWithMeasurement() < minNumberOfLayers) @@ -116,16 +136,16 @@ void SingleLongTrackProducer::produce(edm::Event &iEvent, const edm::EventSetup bestTrack.setExtra(track.extra()); } if (debug) - edm::LogPrint("SingleLongTrackProducer") << " deltaR (general) track to matched Track: " << dRmin << std::endl; + edm::LogPrint("SingleLongTrackProducer") << " deltaR (general) track to matched Track: " << dRmin; if (debug) - edm::LogPrint("SingleLongTrackProducer") << "chi2Ndof:" << chiNdof << " best Track: " << fitProb << std::endl; + edm::LogPrint("SingleLongTrackProducer") << "chi2Ndof:" << chiNdof << " best Track: " << fitProb; } selTracks.push_back(bestTrack); if (debug) - edm::LogPrint("SingleLongTrackProducer") << " number of Tracker Muons: " << tMuon << ", thereof " - << selTracks.size() << " tracks passed preselection." << std::endl; + edm::LogPrint("SingleLongTrackProducer") + << " number of Tracker Muons: " << tMuon << ", thereof " << selTracks.size() << " tracks passed preselection."; // check hits validity in preselected tracks bool hitIsNotValid{false}; @@ -170,14 +190,14 @@ void SingleLongTrackProducer::produce(edm::Event &iEvent, const edm::EventSetup if (onlyValidHits && !hit.isValid()) { if (debug) - edm::LogPrint("SingleLongTrackProducer") << "hit not valid: " << h << std::endl; + edm::LogPrint("SingleLongTrackProducer") << "hit not valid: " << h; continue; } // loop over the hits of the track. if (onlyValidHits && !(hitpattern.validHitFilter(pHit))) { if (debug) - edm::LogPrint("SingleLongTrackProducer") << "hit not valid: " << h << std::endl; + edm::LogPrint("SingleLongTrackProducer") << "hit not valid: " << h; continue; } } @@ -193,6 +213,13 @@ void SingleLongTrackProducer::produce(edm::Event &iEvent, const edm::EventSetup << "found tracks with " << deref << "missing valid hits and " << deref2 << " missing hit pattern"; } + if (debug) { + auto const &moduleType = moduleDescription().moduleName(); + auto const &moduleLabel = moduleDescription().moduleLabel(); + edm::LogPrint("SingleLongTrackProducer") << "[" << moduleType << "] (" << moduleLabel << ") " + << " output track size: " << goodTracks.get()->size(); + } + // save track collection in event: iEvent.put(std::move(goodTracks), ""); } @@ -200,7 +227,7 @@ void SingleLongTrackProducer::produce(edm::Event &iEvent, const edm::EventSetup void SingleLongTrackProducer::fillDescriptions(edm::ConfigurationDescriptions &descriptions) { edm::ParameterSetDescription desc; desc.add("allTracks", edm::InputTag("generalTracks")); - desc.add("matchMuons", edm::InputTag("muons")); + desc.add("matchMuons", edm::InputTag("earlyMuons")); desc.add("PrimaryVertex", edm::InputTag("offlinePrimaryVertices")); desc.add("minNumberOfLayers", 10); desc.add("requiredDr", 0.01); diff --git a/RecoTracker/FinalTrackSelectors/python/SingleLongTrackProducer_cfi.py b/RecoTracker/FinalTrackSelectors/python/SingleLongTrackProducer_cfi.py index 181a536038cfe..1f97050cea9ba 100644 --- a/RecoTracker/FinalTrackSelectors/python/SingleLongTrackProducer_cfi.py +++ b/RecoTracker/FinalTrackSelectors/python/SingleLongTrackProducer_cfi.py @@ -4,7 +4,7 @@ SingleLongTrackProducer = singleLongTrackProducer.clone( allTracks = "generalTracks", - matchMuons = "muons", + matchMuons = "earlyMuons", requiredDr= 0.01, minNumberOfLayers = 10, onlyValidHits = True,