-
Notifications
You must be signed in to change notification settings - Fork 4.6k
Introduce support for cluster splitting and JetCore iteration in Phase-2 track reconstruction #46001
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
Introduce support for cluster splitting and JetCore iteration in Phase-2 track reconstruction #46001
Changes from all commits
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,3 @@ | ||
| import FWCore.ParameterSet.Config as cms | ||
|
|
||
| jetCoreInPhase2 = cms.Modifier() |
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,3 @@ | ||
| import FWCore.ParameterSet.Config as cms | ||
|
|
||
| splitClustersInPhase2Pixel = cms.Modifier() |
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -1,7 +1,9 @@ | ||
| #include "FWCore/Framework/interface/stream/EDProducer.h" | ||
|
|
||
| #include "FWCore/Framework/interface/Event.h" | ||
| #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h" | ||
| #include "FWCore/ParameterSet/interface/ParameterSet.h" | ||
| #include "FWCore/ParameterSet/interface/ParameterSetDescription.h" | ||
| #include "FWCore/Utilities/interface/InputTag.h" | ||
| #include "DataFormats/Common/interface/Handle.h" | ||
| #include "FWCore/Framework/interface/ESHandle.h" | ||
|
|
@@ -12,9 +14,12 @@ | |
| #include "RecoLocalTracker/Records/interface/TkPixelCPERecord.h" | ||
|
|
||
| #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h" | ||
| #include "Geometry/CommonTopologies/interface/PixelTopology.h" | ||
| #include "Geometry/Records/interface/GlobalTrackingGeometryRecord.h" | ||
| #include "DataFormats/GeometryVector/interface/VectorUtil.h" | ||
|
|
||
| #include "DataFormats/TrackerCommon/interface/TrackerTopology.h" | ||
|
|
||
| #include "DataFormats/VertexReco/interface/Vertex.h" | ||
| #include "DataFormats/VertexReco/interface/VertexFwd.h" | ||
| #include "DataFormats/JetReco/interface/Jet.h" | ||
|
|
@@ -25,6 +30,7 @@ | |
|
|
||
| class JetCoreClusterSplitter : public edm::stream::EDProducer<> { | ||
| public: | ||
| static void fillDescriptions(edm::ConfigurationDescriptions& descriptions); | ||
| JetCoreClusterSplitter(const edm::ParameterSet& iConfig); | ||
| ~JetCoreClusterSplitter() override; | ||
| void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override; | ||
|
|
@@ -49,11 +55,17 @@ class JetCoreClusterSplitter : public edm::stream::EDProducer<> { | |
|
|
||
| edm::ESGetToken<GlobalTrackingGeometry, GlobalTrackingGeometryRecord> const tTrackingGeom_; | ||
| edm::ESGetToken<PixelClusterParameterEstimator, TkPixelCPERecord> const tCPE_; | ||
| edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> const tTrackerTopo_; | ||
|
|
||
| bool verbose; | ||
| double ptMin_; | ||
| double deltaR_; | ||
| double chargeFracMin_; | ||
| float expSizeXAtLorentzAngleIncidence_; | ||
| float expSizeXDeltaPerTanAlpha_; | ||
| float expSizeYAtNormalIncidence_; | ||
| float tanLorentzAngle_; | ||
| float tanLorentzAngleBarrelLayer1_; | ||
| edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster>> pixelClusters_; | ||
| edm::EDGetTokenT<reco::VertexCollection> vertices_; | ||
| edm::EDGetTokenT<edm::View<reco::Candidate>> cores_; | ||
|
|
@@ -64,13 +76,45 @@ class JetCoreClusterSplitter : public edm::stream::EDProducer<> { | |
| double centralMIPCharge_; | ||
| }; | ||
|
|
||
| void JetCoreClusterSplitter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) { | ||
| edm::ParameterSetDescription desc; | ||
|
|
||
| desc.add<std::string>("pixelCPE", "PixelCPEGeneric"); | ||
| desc.add<bool>("verbose", false); | ||
| desc.add<double>("ptMin", 200.0); | ||
| desc.add<double>("deltaRmax", 0.05); | ||
| desc.add<double>("chargeFractionMin", 2.0); | ||
| desc.add<edm::InputTag>("pixelClusters", edm::InputTag("siPixelCluster")); | ||
| desc.add<edm::InputTag>("vertices", edm::InputTag("offlinePrimaryVertices")); | ||
| desc.add<edm::InputTag>("cores", edm::InputTag("ak5CaloJets")); | ||
| desc.add<double>("forceXError", 100.0); | ||
| desc.add<double>("forceYError", 150.0); | ||
| desc.add<double>("fractionalWidth", 0.4); | ||
| desc.add<double>("chargePerUnit", 2000.0); | ||
| desc.add<double>("centralMIPCharge", 26e3); | ||
|
|
||
| desc.add<double>("expSizeXAtLorentzAngleIncidence", 1.5); | ||
| desc.add<double>("expSizeXDeltaPerTanAlpha", 0.0); | ||
| desc.add<double>("expSizeYAtNormalIncidence", 1.3); | ||
| desc.add<double>("tanLorentzAngle", 0.0); | ||
| desc.add<double>("tanLorentzAngleBarrelLayer1", 0.0); | ||
|
|
||
| descriptions.addWithDefaultLabel(desc); | ||
| } | ||
|
|
||
| JetCoreClusterSplitter::JetCoreClusterSplitter(const edm::ParameterSet& iConfig) | ||
| : tTrackingGeom_(esConsumes()), | ||
| tCPE_(esConsumes(edm::ESInputTag("", iConfig.getParameter<std::string>("pixelCPE")))), | ||
| tTrackerTopo_(esConsumes()), | ||
| verbose(iConfig.getParameter<bool>("verbose")), | ||
| ptMin_(iConfig.getParameter<double>("ptMin")), | ||
| deltaR_(iConfig.getParameter<double>("deltaRmax")), | ||
| chargeFracMin_(iConfig.getParameter<double>("chargeFractionMin")), | ||
| expSizeXAtLorentzAngleIncidence_(iConfig.getParameter<double>("expSizeXAtLorentzAngleIncidence")), | ||
| expSizeXDeltaPerTanAlpha_(iConfig.getParameter<double>("expSizeXDeltaPerTanAlpha")), | ||
| expSizeYAtNormalIncidence_(iConfig.getParameter<double>("expSizeYAtNormalIncidence")), | ||
| tanLorentzAngle_(iConfig.getParameter<double>("tanLorentzAngle")), | ||
| tanLorentzAngleBarrelLayer1_(iConfig.getParameter<double>("tanLorentzAngleBarrelLayer1")), | ||
| pixelClusters_( | ||
| consumes<edmNew::DetSetVector<SiPixelCluster>>(iConfig.getParameter<edm::InputTag>("pixelClusters"))), | ||
| vertices_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"))), | ||
|
|
@@ -92,6 +136,7 @@ bool SortPixels(const SiPixelCluster::Pixel& i, const SiPixelCluster::Pixel& j) | |
| void JetCoreClusterSplitter::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) { | ||
| using namespace edm; | ||
| const auto& geometry = &iSetup.getData(tTrackingGeom_); | ||
| const auto& topology = &iSetup.getData(tTrackerTopo_); | ||
|
|
||
| Handle<edmNew::DetSetVector<SiPixelCluster>> inputPixelClusters; | ||
| iEvent.getByToken(pixelClusters_, inputPixelClusters); | ||
|
|
@@ -111,6 +156,13 @@ void JetCoreClusterSplitter::produce(edm::Event& iEvent, const edm::EventSetup& | |
| edmNew::DetSetVector<SiPixelCluster>::FastFiller filler(*output, detIt->id()); | ||
| const edmNew::DetSet<SiPixelCluster>& detset = *detIt; | ||
| const GeomDet* det = geometry->idToDet(detset.id()); | ||
| float pitchX, pitchY; | ||
| std::tie(pitchX, pitchY) = static_cast<const PixelTopology&>(det->topology()).pitch(); | ||
| float thickness = det->surface().bounds().thickness(); | ||
| float tanLorentzAngle = tanLorentzAngle_; | ||
| if (DetId(detset.id()).subdetId() == 1 /* px barrel */ && topology->pxbLayer(detset.id()) == 1) { | ||
| tanLorentzAngle = tanLorentzAngleBarrelLayer1_; | ||
|
Contributor
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. shouldn't this come from the database?
Author
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. These parameters were originally hardcoded. I gave them the names they now have to try to interpret them within a simple geometrical model, but I'm not sure how close they are to the quantity stored by the same name in the DB, as for Phase-1 I've tried to preserve to values already in place. As the tracking performance seems to be not too sensitive to the fine-tuning of these parameters (from what I could see), if possible I would keep this implementation as it currently is for this PR, and study the effect of using the parameters from the DB as a future development.
Contributor
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more.
I think it would be nice to study the performance doing a scan of the values around the current parametrization. It might even be that we have sub-par performance with the current Run-3 reconstruction because of that. |
||
| } | ||
| for (auto cluster = detset.begin(); cluster != detset.end(); cluster++) { | ||
| const SiPixelCluster& aCluster = *cluster; | ||
| bool hasBeenSplit = false; | ||
|
|
@@ -126,18 +178,18 @@ void JetCoreClusterSplitter::produce(edm::Event& iEvent, const edm::EventSetup& | |
| if (Geom::deltaR(jetDir, clusterDir) < deltaR_) { | ||
| // check if the cluster has to be splitted | ||
|
|
||
| bool isEndCap = (std::abs(cPos.z()) > 30.f); // FIXME: check detID instead! | ||
| float jetZOverRho = jet.momentum().Z() / jet.momentum().Rho(); | ||
| if (isEndCap) | ||
| jetZOverRho = jet.momentum().Rho() / jet.momentum().Z(); | ||
| float expSizeY = std::sqrt((1.3f * 1.3f) + (1.9f * 1.9f) * jetZOverRho * jetZOverRho); | ||
| LocalVector jetDirLocal = det->surface().toLocal(jetDir); | ||
| float jetTanAlpha = jetDirLocal.x() / jetDirLocal.z(); | ||
| float jetTanBeta = jetDirLocal.y() / jetDirLocal.z(); | ||
| float jetZOverRho = std::sqrt(jetTanAlpha * jetTanAlpha + jetTanBeta * jetTanBeta); | ||
| float expSizeX = expSizeXAtLorentzAngleIncidence_ + | ||
| std::abs(expSizeXDeltaPerTanAlpha_ * (jetTanAlpha - tanLorentzAngle)); | ||
| float expSizeY = std::sqrt((expSizeYAtNormalIncidence_ * expSizeYAtNormalIncidence_) + | ||
| thickness * thickness / (pitchY * pitchY) * jetTanBeta * jetTanBeta); | ||
| if (expSizeX < 1.f) | ||
| expSizeX = 1.f; | ||
| if (expSizeY < 1.f) | ||
| expSizeY = 1.f; | ||
| float expSizeX = 1.5f; | ||
| if (isEndCap) { | ||
|
Contributor
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. why the endcap special treatment is not applied anymore ?
Author
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. In general, the treatment that was applied in the case of a module from the endcaps should now be covered by computing everything in the local frame of the module. However, it is also true that I initially forgot to update the definition of (As a general comment, the value of |
||
| expSizeX = expSizeY; | ||
| expSizeY = 1.5f; | ||
| } // in endcap col/rows are switched | ||
| float expCharge = std::sqrt(1.08f + jetZOverRho * jetZOverRho) * centralMIPCharge_; | ||
|
|
||
| if (aCluster.charge() > expCharge * chargeFracMin_ && | ||
|
|
@@ -163,8 +215,10 @@ void JetCoreClusterSplitter::produce(edm::Event& iEvent, const edm::EventSetup& | |
| if (shouldBeSplit) { | ||
| // blowup the error if we failed to split a splittable cluster (does | ||
| // it ever happen) | ||
| c.setSplitClusterErrorX(c.sizeX() * (100.f / 3.f)); // this is not really blowing up .. TODO: tune | ||
| c.setSplitClusterErrorY(c.sizeY() * (150.f / 3.f)); | ||
| const float fromCentiToMicro = 1e4; | ||
| c.setSplitClusterErrorX(c.sizeX() * | ||
| (pitchX * fromCentiToMicro / 3.f)); // this is not really blowing up .. TODO: tune | ||
| c.setSplitClusterErrorY(c.sizeY() * (pitchY * fromCentiToMicro / 3.f)); | ||
| } | ||
| filler.push_back(c); | ||
| std::push_heap(filler.begin(), filler.end(), [](SiPixelCluster const& cl1, SiPixelCluster const& cl2) { | ||
|
|
||
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
shouldn't this come from the database?