diff --git a/PhysicsTools/PatAlgos/plugins/PATLostTracks.cc b/PhysicsTools/PatAlgos/plugins/PATLostTracks.cc index 6fc798eab4d11..7c5be135a21a2 100644 --- a/PhysicsTools/PatAlgos/plugins/PATLostTracks.cc +++ b/PhysicsTools/PatAlgos/plugins/PATLostTracks.cc @@ -68,6 +68,8 @@ namespace pat { const double maxDxyForNotReconstructedPrimary_; const double maxDxySigForNotReconstructedPrimary_; const bool useLegacySetup_; + const bool xiSelection_; + const double xiMassCut_; }; } // namespace pat @@ -100,7 +102,9 @@ pat::PATLostTracks::PATLostTracks(const edm::ParameterSet& iConfig) .getParameter("maxDxyForNotReconstructedPrimary")), maxDxySigForNotReconstructedPrimary_(iConfig.getParameter("pvAssignment") .getParameter("maxDxySigForNotReconstructedPrimary")), - useLegacySetup_(iConfig.getParameter("useLegacySetup")) { + useLegacySetup_(iConfig.getParameter("useLegacySetup")), + xiSelection_(iConfig.getParameter("xiSelection")), + xiMassCut_(iConfig.getParameter("xiMassCut")) { std::vector trkQuals(iConfig.getParameter>("qualsToAutoAccept")); std::transform( trkQuals.begin(), trkQuals.end(), std::back_inserter(qualsToAutoAccept_), reco::TrackBase::qualityByName); @@ -190,10 +194,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(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 mapping(tracks->size(), -1); diff --git a/PhysicsTools/PatAlgos/python/slimming/lostTracks_cfi.py b/PhysicsTools/PatAlgos/python/slimming/lostTracks_cfi.py index e47fe5bf1db65..93f75eb1cc3f8 100644 --- a/PhysicsTools/PatAlgos/python/slimming/lostTracks_cfi.py +++ b/PhysicsTools/PatAlgos/python/slimming/lostTracks_cfi.py @@ -22,7 +22,9 @@ minPtToStoreLowQualityProps = cms.double(0.0), passThroughCut = cms.string("pt>2"), pvAssignment = primaryVertexAssociation.assignment, - useLegacySetup = cms.bool(False) #When True: check only if track used to fit vertex[0] and do not store track detailed info for Pt between 0.5 and minPtToStoreProps GeV + useLegacySetup = cms.bool(False), #When True: check only if track used to fit vertex[0] and do not store track detailed info for Pt between 0.5 and minPtToStoreProps GeV + xiSelection = cms.bool(True), + xiMassCut = cms.double(1.5) ) from Configuration.Eras.Modifier_phase1Pixel_cff import phase1Pixel phase1Pixel.toModify(lostTracks, covarianceVersion =1 )