Skip to content

Commit cf422a5

Browse files
Added some configurables and changed cuts
1 parent f523055 commit cf422a5

1 file changed

Lines changed: 25 additions & 12 deletions

File tree

PWGLF/Tasks/Resonances/lambda1520Analysispo.cxx

Lines changed: 25 additions & 12 deletions
Original file line numberDiff line numberDiff line change
@@ -94,8 +94,12 @@ struct Lambda1520Analysispo {
9494
// ── Track quality cuts ───────────────────────────────────────────────────
9595
Configurable<float> minTrackPt{
9696
"minTrackPt", 0.15f, "Minimum transverse momentum pT of a track [GeV/c]"};
97-
Configurable<float> minTrackMomentum{
98-
"minTrackMomentum", 0.f, "Minimum total momentum p of a track [GeV/c]"};
97+
Configurable<float> minProtonMomentum{
98+
"minProtonMomentum", 0.f,
99+
"Minimum total momentum p for proton PID [GeV/c]"};
100+
Configurable<float> minKaonMomentum{
101+
"minKaonMomentum", 0.f,
102+
"Minimum total momentum p for kaon PID [GeV/c]"};
99103
Configurable<float> minPseudorapidity{
100104
"minPseudorapidity", -0.8f,
101105
"Minimum pseudorapidity eta (detector acceptance)"};
@@ -363,8 +367,8 @@ struct Lambda1520Analysispo {
363367
Configurable<bool> mcRequireTriggerTVX{
364368
"mcRequireTriggerTVX", false,
365369
"MC event selection: require the TVX (T0 vertex) trigger"};
366-
Configurable<bool> mcRequireRecoINELgt0{
367-
"mcRequireRecoINELgt0", false,
370+
Configurable<bool> cEvtRecINELgt0{
371+
"cEvtRecINELgt0", false,
368372
"MC event selection: require reconstructed INEL>0 condition"};
369373

370374
// ── Histogram registry ───────────────────────────────────────────────────
@@ -769,7 +773,7 @@ struct Lambda1520Analysispo {
769773
if (track.tofNSigmaPr() < minTOFNSigmaProton)
770774
return false;
771775

772-
if (combinedNSigmaCutProton < 0 && totalMomentum >= minTrackMomentum) {
776+
if (combinedNSigmaCutProton < 0 && totalMomentum >= minProtonMomentum) {
773777
for (int i = 0; i < nTOFBins - 1; ++i) {
774778
if (totalMomentum >= tofMomBins[i] &&
775779
totalMomentum < tofMomBins[i + 1] &&
@@ -786,15 +790,15 @@ struct Lambda1520Analysispo {
786790
tpcPIDPassed = true;
787791
}
788792

789-
if ((combinedNSigmaCutProton > 0) && totalMomentum >= minTrackMomentum &&
793+
if ((combinedNSigmaCutProton > 0) && totalMomentum >= minProtonMomentum &&
790794
(combinedNSigProton < circularCutSquared &&
791795
combinedNSigPion > circularRejCutSquared &&
792796
combinedNSigKaon > circularRejCutSquared)) {
793797
tofPIDPassed = true;
794798
tpcPIDPassed = true;
795799
}
796800

797-
if (totalMomentum < minTrackMomentum &&
801+
if (totalMomentum < minProtonMomentum &&
798802
tpcNSigProton < maxTPCNSigmaProton) {
799803
tofPIDPassed = true;
800804
tpcPIDPassed = true;
@@ -853,7 +857,7 @@ struct Lambda1520Analysispo {
853857
if (track.tofNSigmaKa() < minTOFNSigmaKaon)
854858
return false;
855859

856-
if (combinedNSigmaCutKaon < 0 && totalMomentum >= minTrackMomentum) {
860+
if (combinedNSigmaCutKaon < 0 && totalMomentum >= minKaonMomentum) {
857861
for (int i = 0; i < nTOFBins - 1; ++i) {
858862
if (totalMomentum >= tofMomBins[i] &&
859863
totalMomentum < tofMomBins[i + 1] &&
@@ -870,15 +874,15 @@ struct Lambda1520Analysispo {
870874
tpcPIDPassed = true;
871875
}
872876

873-
if ((combinedNSigmaCutKaon > 0) && totalMomentum >= minTrackMomentum &&
877+
if ((combinedNSigmaCutKaon > 0) && totalMomentum >= minKaonMomentum &&
874878
(combinedNSigKaon < circularCutSquared &&
875879
combinedNSigPion > circularRejCutSquared &&
876880
combinedNSigProton > circularRejCutSquared)) {
877881
tofPIDPassed = true;
878882
tpcPIDPassed = true;
879883
}
880884

881-
if (totalMomentum < minTrackMomentum && tpcNSigKaon < maxTPCNSigmaKaon) {
885+
if (totalMomentum < minKaonMomentum && tpcNSigKaon < maxTPCNSigmaKaon) {
882886
tofPIDPassed = true;
883887
tpcPIDPassed = true;
884888
}
@@ -1236,6 +1240,9 @@ struct Lambda1520Analysispo {
12361240
void processData(ResonanceCollisionsWithEP::iterator const& collision,
12371241
ResonanceTrackTable const& tracks)
12381242
{
1243+
if (cEvtRecINELgt0 && !collision.isRecINELgt0()) // Check reco INELgt0 (at least one PV track in |eta| < 1) about the collision
1244+
return;
1245+
12391246
allHistograms.fill(HIST("Event/centralityVsOccupancy"), collision.cent(),
12401247
100);
12411248
allHistograms.fill(HIST("Event/primaryVertexZ"), collision.posZ());
@@ -1270,7 +1277,7 @@ struct Lambda1520Analysispo {
12701277
return;
12711278
allHistograms.fill(HIST("Event/mcEventSelectionCutflow"), 4);
12721279

1273-
if (mcRequireRecoINELgt0 && !collision.isRecINELgt0())
1280+
if (cEvtRecINELgt0 && !collision.isRecINELgt0())
12741281
return;
12751282
allHistograms.fill(HIST("Event/mcEventSelectionCutflow"), 5);
12761283

@@ -1370,7 +1377,7 @@ struct Lambda1520Analysispo {
13701377
if (mcRequireSel8 && !collision.isInSel8())
13711378
return;
13721379
allHistograms.fill(HIST("SignalLoss/mcEventSelectionCutflow"), 4);
1373-
if (mcRequireRecoINELgt0 && !collision.isRecINELgt0())
1380+
if (cEvtRecINELgt0 && !collision.isRecINELgt0())
13741381
return;
13751382
allHistograms.fill(HIST("SignalLoss/mcEventSelectionCutflow"), 5);
13761383
if (mcRequireAfterAllCuts && !collision.isInAfterAllCuts())
@@ -1479,6 +1486,8 @@ struct Lambda1520Analysispo {
14791486

14801487
// FIX const-ref-in-for-loop: const auto& structured binding
14811488
for (const auto& [col1, tracks1, col2, tracks2] : eventPairs) {
1489+
if (cEvtRecINELgt0 && !col1.isRecINELgt0()) // Check reco INELgt0 (at least one PV track in |eta| < 1) about the collision
1490+
return;
14821491
allHistograms.fill(HIST("Event/mixingBins_centralityVsVtxZVsEventPlane"),
14831492
col1.cent(), col1.posZ(), col1.evtPl());
14841493
fillInvariantMassHistograms<true, false>(tracks1, tracks2, col1.cent());
@@ -1503,6 +1512,8 @@ struct Lambda1520Analysispo {
15031512
if (doprocessData)
15041513
LOG(error) << "Disable processData() first when using processDatadf()!";
15051514

1515+
if (cEvtRecINELgt0 && !collision.isRecINELgt0()) // Check reco INELgt0 (at least one PV track in |eta| < 1) about the collision
1516+
return;
15061517
auto occupancyValue = 100;
15071518
if (applyOccupancyInTimeRangeCut)
15081519
occupancyValue = collision.trackOccupancyInTimeRange();
@@ -1538,6 +1549,8 @@ struct Lambda1520Analysispo {
15381549

15391550
// FIX const-ref-in-for-loop: const auto& structured binding
15401551
for (const auto& [col1, tracks1, col2, tracks2] : eventPairs) {
1552+
if (cEvtRecINELgt0 && !col1.isRecINELgt0()) // Check reco INELgt0 (at least one PV track in |eta| < 1) about the collision
1553+
return;
15411554
auto occupancyValue = 100;
15421555
if (applyOccupancyInTimeRangeCut)
15431556
occupancyValue = col1.trackOccupancyInTimeRange();

0 commit comments

Comments
 (0)