Skip to content

Commit ac96e41

Browse files
committed
Fix o2linter warnings on TLorenzVector
1 parent e42faca commit ac96e41

1 file changed

Lines changed: 36 additions & 40 deletions

File tree

PWGUD/Tasks/FwdMuonsUPC.cxx

Lines changed: 36 additions & 40 deletions
Original file line numberDiff line numberDiff line change
@@ -23,8 +23,9 @@
2323
#include "Framework/O2DatabasePDGPlugin.h"
2424
#include "Framework/runDataProcessing.h"
2525

26-
#include "TLorentzVector.h"
2726
#include "TRandom3.h"
27+
#include <Math/Vector4D.h>
28+
#include <Math/VectorUtil.h>
2829

2930
#include <unordered_map>
3031
#include <vector>
@@ -349,9 +350,9 @@ struct FwdMuonsUPC {
349350
{
350351
float rAbs = fwdTrack.rAtAbsorberEnd();
351352
float pDca = fwdTrack.pDca();
352-
TLorentzVector p;
353353
auto mMu = particleMass(kMuonPDG);
354-
p.SetXYZM(fwdTrack.px(), fwdTrack.py(), fwdTrack.pz(), mMu);
354+
ROOT::Math::PxPyPzMVector p{fwdTrack.px(), fwdTrack.py(), fwdTrack.pz(), mMu};
355+
355356
float eta = p.Eta();
356357
float pt = p.Pt();
357358
float pDcaMax = rAbs < kRAbsMid ? kPDca1 : kPDca2;
@@ -368,9 +369,9 @@ struct FwdMuonsUPC {
368369
}
369370

370371
// function to compute phi for azimuth anisotropy
371-
void computePhiAnis(TLorentzVector p1, TLorentzVector p2, int sign1, float& phiAverage, float& phiCharge)
372+
void computePhiAnis(ROOT::Math::PxPyPzMVector p1, ROOT::Math::PxPyPzMVector p2, int sign1, float& phiAverage, float& phiCharge)
372373
{
373-
TLorentzVector tSum, tDiffAv, tDiffCh;
374+
ROOT::Math::PxPyPzMVector tSum, tDiffAv, tDiffCh;
374375
tSum = p1 + p2;
375376
float halfUnity = 0.5;
376377
if (sign1 > 0) {
@@ -388,9 +389,9 @@ struct FwdMuonsUPC {
388389
}
389390

390391
// average
391-
phiAverage = tSum.DeltaPhi(tDiffAv);
392+
phiAverage = ROOT::Math::VectorUtil::DeltaPhi(tSum, tDiffAv);
392393
// charge
393-
phiCharge = tSum.DeltaPhi(tDiffCh);
394+
phiCharge = ROOT::Math::VectorUtil::DeltaPhi(tSum, tDiffCh);
394395
}
395396

396397
// function that processes the candidates:
@@ -453,11 +454,10 @@ struct FwdMuonsUPC {
453454
return;
454455

455456
// form Lorentz vectors
456-
TLorentzVector p1, p2;
457457
auto mMu = particleMass(kMuonPDG);
458-
p1.SetXYZM(tr1.px(), tr1.py(), tr1.pz(), mMu);
459-
p2.SetXYZM(tr2.px(), tr2.py(), tr2.pz(), mMu);
460-
TLorentzVector p = p1 + p2;
458+
ROOT::Math::PxPyPzMVector p1{tr1.px(), tr1.py(), tr1.pz(), mMu};
459+
ROOT::Math::PxPyPzMVector p2{tr2.px(), tr2.py(), tr2.pz(), mMu};
460+
ROOT::Math::PxPyPzMVector p = p1 + p2;
461461

462462
// cut on pair kinematics
463463
// select mass
@@ -520,15 +520,15 @@ struct FwdMuonsUPC {
520520
dimuSel(cand.runNumber(),
521521
p.M(), p.E(), p.Px(), p.Py(), p.Pz(), p.Pt(), p.Rapidity(), p.Phi(),
522522
phiAverage, phiCharge,
523-
p1.E(), p1.Px(), p1.Py(), p1.Pz(), p1.Pt(), p1.PseudoRapidity(), p1.Phi(), static_cast<int>(myTrackType),
524-
p2.E(), p2.Px(), p2.Py(), p2.Pz(), p2.Pt(), p2.PseudoRapidity(), p2.Phi(), static_cast<int>(myTrackType),
523+
p1.E(), p1.Px(), p1.Py(), p1.Pz(), p1.Pt(), p1.Eta(), p1.Phi(), static_cast<int>(myTrackType),
524+
p2.E(), p2.Px(), p2.Py(), p2.Pz(), p2.Pt(), p2.Eta(), p2.Phi(), static_cast<int>(myTrackType),
525525
zdc.timeA, zdc.enA, zdc.timeC, zdc.enC, znClass);
526526
} else {
527527
dimuSel(cand.runNumber(),
528528
p.M(), p.E(), p.Px(), p.Py(), p.Pz(), p.Pt(), p.Rapidity(), p.Phi(),
529529
phiAverage, phiCharge,
530-
p2.E(), p2.Px(), p2.Py(), p2.Pz(), p2.Pt(), p2.PseudoRapidity(), p2.Phi(), static_cast<int>(myTrackType),
531-
p1.E(), p1.Px(), p1.Py(), p1.Pz(), p1.Pt(), p1.PseudoRapidity(), p1.Phi(), static_cast<int>(myTrackType),
530+
p2.E(), p2.Px(), p2.Py(), p2.Pz(), p2.Pt(), p2.Eta(), p2.Phi(), static_cast<int>(myTrackType),
531+
p1.E(), p1.Px(), p1.Py(), p1.Pz(), p1.Pt(), p1.Eta(), p1.Phi(), static_cast<int>(myTrackType),
532532
zdc.timeA, zdc.enA, zdc.timeC, zdc.enC, znClass);
533533
}
534534
}
@@ -549,11 +549,10 @@ struct FwdMuonsUPC {
549549
}
550550

551551
// create Lorentz vectors
552-
TLorentzVector p1, p2;
553552
auto mMu = particleMass(kMuonPDG);
554-
p1.SetXYZM(McPart1.px(), McPart1.py(), McPart1.pz(), mMu);
555-
p2.SetXYZM(McPart2.px(), McPart2.py(), McPart2.pz(), mMu);
556-
TLorentzVector p = p1 + p2;
553+
ROOT::Math::PxPyPzMVector p1{McPart1.px(), McPart1.py(), McPart1.pz(), mMu};
554+
ROOT::Math::PxPyPzMVector p2{McPart2.px(), McPart2.py(), McPart2.pz(), mMu};
555+
ROOT::Math::PxPyPzMVector p = p1 + p2;
557556

558557
// cut on pair kinematics
559558
// select mass
@@ -587,13 +586,13 @@ struct FwdMuonsUPC {
587586
if (McPart1.pdgCode() < 0) {
588587
dimuGen(p.M(), p.Pt(), p.Rapidity(), p.Phi(),
589588
phiAverage, phiCharge,
590-
p1.Pt(), p1.PseudoRapidity(), p1.Phi(),
591-
p2.Pt(), p2.PseudoRapidity(), p2.Phi());
589+
p1.Pt(), p1.Eta(), p1.Phi(),
590+
p2.Pt(), p2.Eta(), p2.Phi());
592591
} else {
593592
dimuGen(p.M(), p.Pt(), p.Rapidity(), p.Phi(),
594593
phiAverage, phiCharge,
595-
p2.Pt(), p2.PseudoRapidity(), p2.Phi(),
596-
p1.Pt(), p1.PseudoRapidity(), p1.Phi());
594+
p2.Pt(), p2.Eta(), p2.Phi(),
595+
p1.Pt(), p1.Eta(), p1.Phi());
597596
}
598597
}
599598

@@ -655,11 +654,10 @@ struct FwdMuonsUPC {
655654
return;
656655

657656
// form Lorentz vectors
658-
TLorentzVector p1, p2;
659657
auto mMu = particleMass(kMuonPDG);
660-
p1.SetXYZM(tr1.px(), tr1.py(), tr1.pz(), mMu);
661-
p2.SetXYZM(tr2.px(), tr2.py(), tr2.pz(), mMu);
662-
TLorentzVector p = p1 + p2;
658+
ROOT::Math::PxPyPzMVector p1{tr1.px(), tr1.py(), tr1.pz(), mMu};
659+
ROOT::Math::PxPyPzMVector p2{tr2.px(), tr2.py(), tr2.pz(), mMu};
660+
ROOT::Math::PxPyPzMVector p = p1 + p2;
663661

664662
// cut on pair kinematics (reco candidates)
665663
// select mass
@@ -683,11 +681,9 @@ struct FwdMuonsUPC {
683681
float phiCharge = 0;
684682
computePhiAnis(p1, p2, tr1.sign(), phiAverage, phiCharge);
685683

686-
// gen particle
687-
TLorentzVector p1Mc, p2Mc;
688-
p1Mc.SetXYZM(McPart1.px(), McPart1.py(), McPart1.pz(), mMu);
689-
p2Mc.SetXYZM(McPart2.px(), McPart2.py(), McPart2.pz(), mMu);
690-
TLorentzVector pMc = p1Mc + p2Mc;
684+
ROOT::Math::PxPyPzMVector p1Mc{McPart1.px(), McPart1.py(), McPart1.pz(), mMu};
685+
ROOT::Math::PxPyPzMVector p2Mc{McPart2.px(), McPart2.py(), McPart2.pz(), mMu};
686+
ROOT::Math::PxPyPzMVector pMc = p2Mc + p2Mc;
691687

692688
// compute gen phi for azimuth anisotropy
693689
float phiGenAverage = 0;
@@ -713,22 +709,22 @@ struct FwdMuonsUPC {
713709
dimuReco(cand.runNumber(),
714710
p.M(), p.Pt(), p.Rapidity(), p.Phi(),
715711
phiAverage, phiCharge,
716-
p1.Pt(), p1.PseudoRapidity(), p1.Phi(), static_cast<int>(myTrackType),
717-
p2.Pt(), p2.PseudoRapidity(), p2.Phi(), static_cast<int>(myTrackType),
712+
p1.Pt(), p1.Eta(), p1.Phi(), static_cast<int>(myTrackType),
713+
p2.Pt(), p2.Eta(), p2.Phi(), static_cast<int>(myTrackType),
718714
// gen info
719715
pMc.Pt(), pMc.Rapidity(), pMc.Phi(),
720-
p1Mc.Pt(), p1Mc.PseudoRapidity(), p1Mc.Phi(),
721-
p2Mc.Pt(), p2Mc.PseudoRapidity(), p2Mc.Phi());
716+
p1Mc.Pt(), p1Mc.Eta(), p1Mc.Phi(),
717+
p2Mc.Pt(), p2Mc.Eta(), p2Mc.Phi());
722718
} else {
723719
dimuReco(cand.runNumber(),
724720
p.M(), p.Pt(), p.Rapidity(), p.Phi(),
725721
phiAverage, phiCharge,
726-
p2.Pt(), p2.PseudoRapidity(), p2.Phi(), static_cast<int>(myTrackType),
727-
p1.Pt(), p1.PseudoRapidity(), p1.Phi(), static_cast<int>(myTrackType),
722+
p2.Pt(), p2.Eta(), p2.Phi(), static_cast<int>(myTrackType),
723+
p1.Pt(), p1.Eta(), p1.Phi(), static_cast<int>(myTrackType),
728724
// gen info
729725
pMc.Pt(), pMc.Rapidity(), pMc.Phi(),
730-
p2Mc.Pt(), p2Mc.PseudoRapidity(), p2Mc.Phi(),
731-
p1Mc.Pt(), p1Mc.PseudoRapidity(), p1Mc.Phi());
726+
p2Mc.Pt(), p2Mc.Eta(), p2Mc.Phi(),
727+
p1Mc.Pt(), p1Mc.Eta(), p1Mc.Phi());
732728
}
733729
}
734730

0 commit comments

Comments
 (0)