Skip to content

Commit 4f044f8

Browse files
author
AdrianoDee
committed
Moving to XYZTLorentz vectors
1 parent 40e810e commit 4f044f8

File tree

2 files changed

+6
-6
lines changed

2 files changed

+6
-6
lines changed

PhysicsTools/PatAlgos/plugins/PATLostTracks.cc

Lines changed: 5 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -105,7 +105,7 @@ pat::PATLostTracks::PATLostTracks(const edm::ParameterSet& iConfig)
105105
maxDxySigForNotReconstructedPrimary_(iConfig.getParameter<edm::ParameterSet>("pvAssignment")
106106
.getParameter<double>("maxDxySigForNotReconstructedPrimary")),
107107
useLegacySetup_(iConfig.getParameter<bool>("useLegacySetup")),
108-
xiSelection_(iConfig.getParameter<double>("xiSelection")),
108+
xiSelection_(iConfig.getParameter<bool>("xiSelection")),
109109
xiMassCut_(iConfig.getParameter<double>("xiMassCut")) {
110110
std::vector<std::string> trkQuals(iConfig.getParameter<std::vector<std::string>>("qualsToAutoAccept"));
111111
std::transform(
@@ -205,14 +205,14 @@ void pat::PATLostTracks::produce(edm::StreamID, edm::Event& iEvent, const edm::E
205205
}
206206
if (xiSelection_) {
207207
// selecting potential Xi- -> Lambda pi candidates
208-
TLorentzVector p4Lambda;
209-
p4Lambda.SetPtEtaPhiM(v0.pt(), v0.eta(), v0.phi(), v0.mass());
208+
math::XYZTLorentzVector p4Lambda(v0.px(), v0.py(), v0.pz(), sqrt(v0.momentum().mag2() + v0.mass()*v0.mass()));
209+
210210
for (unsigned int trkIndx = 0; trkIndx < tracks->size(); trkIndx++) {
211211
reco::TrackRef trk(tracks, trkIndx);
212212
if ((*trk).charge() * protonCharge < 0)
213213
continue;
214-
TLorentzVector p4pi;
215-
p4pi.SetPtEtaPhiM((*trk).pt(), (*trk).eta(), (*trk).phi(), 0.13957061);
214+
215+
math::XYZTLorentzVector p4pi(trk->px(), trk->py(), trk->pz(), trk->momentum().mag2() + 0.01947995518); // pion mass ^2
216216
if ((p4Lambda + p4pi).M() < xiMassCut_) { // selecting potential Xi- candidates
217217
if (trkStatus[trkIndx] == TrkStatus::NOTUSED)
218218
trkStatus[trkIndx] = TrkStatus::VTX;

PhysicsTools/PatAlgos/python/slimming/lostTracks_cfi.py

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -24,7 +24,7 @@
2424
pvAssignment = primaryVertexAssociation.assignment,
2525
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
2626
xiSelection = cms.bool(True),
27-
xiMassCut = cms.bool(1.5)
27+
xiMassCut = cms.double(1.5)
2828
)
2929
from Configuration.Eras.Modifier_phase1Pixel_cff import phase1Pixel
3030
phase1Pixel.toModify(lostTracks, covarianceVersion =1 )

0 commit comments

Comments
 (0)