Skip to content

Commit 6ead9cd

Browse files
authored
[PWGJE] add min pThat cut (AliceO2Group#9603)
1 parent c21bef5 commit 6ead9cd

File tree

3 files changed

+24
-9
lines changed

3 files changed

+24
-9
lines changed

PWGJE/Tasks/jetFinderFullQA.cxx

Lines changed: 8 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -74,6 +74,7 @@ struct JetFinderFullQATask {
7474
Configurable<float> pTHatMaxMCD{"pTHatMaxMCD", 999.0, "maximum fraction of hard scattering for jet acceptance in detector MC"};
7575
Configurable<float> pTHatMaxMCP{"pTHatMaxMCP", 999.0, "maximum fraction of hard scattering for jet acceptance in particle MC"};
7676
Configurable<float> pTHatExponent{"pTHatExponent", 6.0, "exponent of the event weight for the calculation of pTHat"};
77+
Configurable<float> pTHatAbsoluteMin{"pTHatAbsoluteMin", -99.0, "minimum value of pTHat"};
7778
Configurable<bool> skipMBGapEvents{"skipMBGapEvents", false, "flag to choose to reject min. bias gap events; jet-level rejection applied at the jet finder level, here rejection is applied for collision and track process functions"};
7879

7980
std::vector<bool> filledJetR;
@@ -278,7 +279,7 @@ struct JetFinderFullQATask {
278279
{
279280

280281
float pTHat = 10. / (std::pow(weight, 1.0 / pTHatExponent));
281-
if (jet.pt() > pTHatMaxMCD * pTHat) {
282+
if (jet.pt() > pTHatMaxMCD * pTHat || pTHat < pTHatAbsoluteMin) {
282283
return;
283284
}
284285
registry.fill(HIST("h_jet_phat_weighted"), pTHat, weight);
@@ -328,7 +329,7 @@ struct JetFinderFullQATask {
328329
{
329330

330331
float pTHat = 10. / (std::pow(weight, 1.0 / pTHatExponent));
331-
if (jet.pt() > pTHatMaxMCP * pTHat) {
332+
if (jet.pt() > pTHatMaxMCP * pTHat || pTHat < pTHatAbsoluteMin) {
332333
return;
333334
}
334335
registry.fill(HIST("h_jet_phat_part_weighted"), pTHat, weight);
@@ -363,7 +364,7 @@ struct JetFinderFullQATask {
363364
{
364365

365366
float pTHat = 10. / (std::pow(weight, 1.0 / pTHatExponent));
366-
if (jetBase.pt() > pTHatMaxMCD * pTHat) {
367+
if (jetBase.pt() > pTHatMaxMCD * pTHat || pTHat < pTHatAbsoluteMin) {
367368
return;
368369
}
369370

@@ -436,6 +437,10 @@ struct JetFinderFullQATask {
436437
template <typename T, typename U>
437438
void fillTrackHistograms(T const& tracks, U const& clusters, float weight = 1.0)
438439
{
440+
float pTHat = 10. / (std::pow(weight, 1.0 / pTHatExponent));
441+
if (pTHat < pTHatAbsoluteMin) {
442+
return;
443+
}
439444
for (auto const& track : tracks) {
440445
if (!jetderiveddatautilities::selectTrack(track, trackSelection)) {
441446
continue;

PWGJE/Tasks/jetFinderHFQA.cxx

Lines changed: 8 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -68,6 +68,7 @@ struct JetFinderHFQATask {
6868
Configurable<float> pTHatMaxMCD{"pTHatMaxMCD", 999.0, "maximum fraction of hard scattering for jet acceptance in detector MC"};
6969
Configurable<float> pTHatMaxMCP{"pTHatMaxMCP", 999.0, "maximum fraction of hard scattering for jet acceptance in particle MC"};
7070
Configurable<float> pTHatExponent{"pTHatExponent", 6.0, "exponent of the event weight for the calculation of pTHat"};
71+
Configurable<float> pTHatAbsoluteMin{"pTHatAbsoluteMin", -99.0, "minimum value of pTHat"};
7172
Configurable<float> randomConeR{"randomConeR", 0.4, "size of random Cone for estimating background fluctuations"};
7273
Configurable<bool> skipMBGapEvents{"skipMBGapEvents", false, "flag to choose to reject min. bias gap events; jet-level rejection applied at the jet finder level, here rejection is applied for collision and track process functions"};
7374

@@ -527,7 +528,7 @@ struct JetFinderHFQATask {
527528
{
528529

529530
float pTHat = 10. / (std::pow(weight, 1.0 / pTHatExponent));
530-
if (jet.pt() > pTHatMaxMCD * pTHat) {
531+
if (jet.pt() > pTHatMaxMCD * pTHat || pTHat < pTHatAbsoluteMin) {
531532
return;
532533
}
533534
registry.fill(HIST("h_jet_phat_weighted"), pTHat, weight);
@@ -656,7 +657,7 @@ struct JetFinderHFQATask {
656657
{
657658

658659
float pTHat = 10. / (std::pow(weight, 1.0 / pTHatExponent));
659-
if (jet.pt() > pTHatMaxMCP * pTHat) {
660+
if (jet.pt() > pTHatMaxMCP * pTHat || pTHat < pTHatAbsoluteMin) {
660661
return;
661662
}
662663
registry.fill(HIST("h_jet_phat_part_weighted"), pTHat, weight);
@@ -698,7 +699,7 @@ struct JetFinderHFQATask {
698699
{
699700

700701
float pTHat = 10. / (std::pow(weight, 1.0 / pTHatExponent));
701-
if (jetBase.pt() > pTHatMaxMCD * pTHat) {
702+
if (jetBase.pt() > pTHatMaxMCD * pTHat || pTHat < pTHatAbsoluteMin) {
702703
return;
703704
}
704705
auto candidateBasePt = 0.0;
@@ -930,6 +931,10 @@ struct JetFinderHFQATask {
930931
template <typename T, typename U>
931932
void fillTrackHistograms(T const& collision, U const& tracks, float weight = 1.0)
932933
{
934+
float pTHat = 10. / (std::pow(weight, 1.0 / pTHatExponent));
935+
if (pTHat < pTHatAbsoluteMin) {
936+
return;
937+
}
933938
for (auto const& track : tracks) {
934939
if (!jetderiveddatautilities::selectTrack(track, trackSelection)) {
935940
continue;

PWGJE/Tasks/jetFinderQA.cxx

Lines changed: 8 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -60,6 +60,7 @@ struct JetFinderQATask {
6060
Configurable<float> pTHatMaxMCD{"pTHatMaxMCD", 999.0, "maximum fraction of hard scattering for jet acceptance in detector MC"};
6161
Configurable<float> pTHatMaxMCP{"pTHatMaxMCP", 999.0, "maximum fraction of hard scattering for jet acceptance in particle MC"};
6262
Configurable<float> pTHatExponent{"pTHatExponent", 6.0, "exponent of the event weight for the calculation of pTHat"};
63+
Configurable<float> pTHatAbsoluteMin{"pTHatAbsoluteMin", -99.0, "minimum value of pTHat"};
6364
Configurable<double> jetPtMax{"jetPtMax", 200., "set jet pT bin max"};
6465
Configurable<float> jetEtaMin{"jetEtaMin", -99.0, "minimum jet pseudorapidity"};
6566
Configurable<float> jetEtaMax{"jetEtaMax", 99.0, "maximum jet pseudorapidity"};
@@ -382,7 +383,7 @@ struct JetFinderQATask {
382383
{
383384

384385
float pTHat = 10. / (std::pow(weight, 1.0 / pTHatExponent));
385-
if (jet.pt() > pTHatMaxMCD * pTHat) {
386+
if (jet.pt() > pTHatMaxMCD * pTHat || pTHat < pTHatAbsoluteMin) {
386387
return;
387388
}
388389
registry.fill(HIST("h_jet_phat"), pTHat);
@@ -482,7 +483,7 @@ struct JetFinderQATask {
482483
{
483484

484485
float pTHat = 10. / (std::pow(weight, 1.0 / pTHatExponent));
485-
if (jet.pt() > pTHatMaxMCP * pTHat) {
486+
if (jet.pt() > pTHatMaxMCP * pTHat || pTHat < pTHatAbsoluteMin) {
486487
return;
487488
}
488489
registry.fill(HIST("h_jet_phat_part"), pTHat);
@@ -512,7 +513,7 @@ struct JetFinderQATask {
512513
void fillMatchedHistograms(T const& jetBase, float weight = 1.0)
513514
{
514515
float pTHat = 10. / (std::pow(weight, 1.0 / pTHatExponent));
515-
if (jetBase.pt() > pTHatMaxMCD * pTHat) {
516+
if (jetBase.pt() > pTHatMaxMCD * pTHat || pTHat < pTHatAbsoluteMin) {
516517
return;
517518
}
518519

@@ -592,6 +593,10 @@ struct JetFinderQATask {
592593
template <typename T, typename U>
593594
void fillTrackHistograms(T const& collision, U const& tracks, float weight = 1.0)
594595
{
596+
float pTHat = 10. / (std::pow(weight, 1.0 / pTHatExponent));
597+
if (pTHat < pTHatAbsoluteMin) {
598+
return;
599+
}
595600
for (auto const& track : tracks) {
596601
if (!jetderiveddatautilities::selectTrack(track, trackSelection)) {
597602
continue;

0 commit comments

Comments
 (0)