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

[Backport] Adding Xi-tracks to lostTracks to 10_6_X #33962

Merged
merged 2 commits into from Jun 4, 2021
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
23 changes: 22 additions & 1 deletion PhysicsTools/PatAlgos/plugins/PATLostTracks.cc
Expand Up @@ -70,6 +70,8 @@ namespace pat {
const double maxDxySigForNotReconstructedPrimary_;
const bool useLegacySetup_;
const bool fillLostInnerHits_;
const bool xiSelection_;
const double xiMassCut_;
};
} // namespace pat

Expand Down Expand Up @@ -104,7 +106,9 @@ pat::PATLostTracks::PATLostTracks(const edm::ParameterSet& iConfig)
maxDxySigForNotReconstructedPrimary_(iConfig.getParameter<edm::ParameterSet>("pvAssignment")
.getParameter<double>("maxDxySigForNotReconstructedPrimary")),
useLegacySetup_(iConfig.getParameter<bool>("useLegacySetup")),
fillLostInnerHits_(iConfig.getParameter<bool>("fillLostInnerHits")) {
fillLostInnerHits_(iConfig.getParameter<bool>("fillLostInnerHits")),
xiSelection_(iConfig.getParameter<bool>("xiSelection")),
xiMassCut_(iConfig.getParameter<double>("xiMassCut")) {
std::vector<std::string> trkQuals(iConfig.getParameter<std::vector<std::string>>("qualsToAutoAccept"));
std::transform(
trkQuals.begin(), trkQuals.end(), std::back_inserter(qualsToAutoAccept_), reco::TrackBase::qualityByName);
Expand Down Expand Up @@ -194,10 +198,27 @@ void pat::PATLostTracks::produce(edm::StreamID, edm::Event& iEvent, const edm::E
}
}
for (const auto& v0 : *lambdas) {
double protonCharge = 0;
for (size_t dIdx = 0; dIdx < v0.numberOfDaughters(); dIdx++) {
size_t key = (dynamic_cast<const reco::RecoChargedCandidate*>(v0.daughter(dIdx)))->track().key();
if (trkStatus[key] == TrkStatus::NOTUSED)
trkStatus[key] = TrkStatus::VTX;
protonCharge += v0.daughter(dIdx)->charge() * v0.daughter(dIdx)->momentum().mag2();
}
if (xiSelection_) {
// selecting potential Xi- -> Lambda pi candidates
math::XYZTLorentzVector p4Lambda(v0.px(), v0.py(), v0.pz(), sqrt(v0.momentum().mag2() + v0.mass() * v0.mass()));

for (unsigned int trkIndx = 0; trkIndx < tracks->size(); trkIndx++) {
reco::TrackRef trk(tracks, trkIndx);
if ((*trk).charge() * protonCharge < 0 || trkStatus[trkIndx] != TrkStatus::NOTUSED)
continue;

math::XYZTLorentzVector p4pi(
trk->px(), trk->py(), trk->pz(), sqrt(trk->momentum().mag2() + 0.01947995518)); // pion mass ^2
if ((p4Lambda + p4pi).M() < xiMassCut_)
trkStatus[trkIndx] = TrkStatus::VTX; // selecting potential Xi- candidates
}
}
}
std::vector<int> mapping(tracks->size(), -1);
Expand Down
6 changes: 4 additions & 2 deletions PhysicsTools/PatAlgos/python/slimming/lostTracks_cfi.py
Expand Up @@ -24,7 +24,9 @@
allowMuonId = cms.bool(False),
pvAssignment = primaryVertexAssociation.assignment,
useLegacySetup = cms.bool(True),
fillLostInnerHits = cms.bool(False)
fillLostInnerHits = cms.bool(False),
xiSelection = cms.bool(False),
xiMassCut = cms.double(1.5)
)
from Configuration.Eras.Modifier_phase1Pixel_cff import phase1Pixel
phase1Pixel.toModify(lostTracks, covarianceVersion =1 )
Expand All @@ -36,4 +38,4 @@
bParking.toModify(lostTracks, fillLostInnerHits = True)

from Configuration.Eras.Modifier_run2_miniAOD_devel_cff import run2_miniAOD_devel
(bParking | run2_miniAOD_devel).toModify(lostTracks, minPtToStoreLowQualityProps = 0.0)
(bParking | run2_miniAOD_devel).toModify(lostTracks, minPtToStoreLowQualityProps = 0.0, xiSelection = True)