Skip to content

Commit 8a4a641

Browse files
authored
[PWGHF] Protect daughter index access in selectEvent to prevent segmentation faults (#2281)
1 parent 92074d3 commit 8a4a641

File tree

1 file changed

+34
-6
lines changed

1 file changed

+34
-6
lines changed

MC/config/PWGHF/external/generator/generator_pythia8_hfhadron_to_nuclei.C

Lines changed: 34 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -3,6 +3,7 @@
33
#include "Pythia8/Pythia.h"
44
#include "TRandom.h"
55
#include <fairlogger/Logger.h>
6+
#include <algorithm>
67
#include <string>
78
#include <vector>
89

@@ -119,18 +120,45 @@ class GeneratorPythia8HFHadToNuclei : public o2::eventgen::GeneratorPythia8
119120
int id = std::abs(event[iPart].id());
120121
float rap = event[iPart].y();
121122
if (id == mHadronPdg && rap > mHadRapidityMin && rap < mHadRapidityMax) {
123+
124+
int d1 = event[iPart].daughter1();
125+
int d2 = event[iPart].daughter2();
126+
127+
// Protection for invalid index of daughter 1
128+
if (d1 < 0 || d1 >= event.size()) {
129+
LOG(info) << "Invalid daughter index 1! d1 = " << d1;
130+
continue;
131+
}
132+
133+
// Protection for invalid index of daughter 2
134+
if (d2 < 0 || d2 >= event.size()) {
135+
LOG(info) << "Invalid daughter index 2! d2 = " << d2;
136+
continue;
137+
}
138+
139+
// Swap order
140+
if (d1 > d2) {
141+
std::swap(d1, d2);
142+
}
143+
144+
// No daughters
145+
if (d1 == 0 && d2 == 0) {
146+
LOG(info) << "No daughters (0,0). Skipping.";
147+
continue;
148+
}
149+
122150
LOG(debug) << "-----------------------------------------------------";
123-
LOG(debug) << "Found hadron " << event[iPart].id() << " with rapidity " << rap << " and daughters " << event[iPart].daughter1() << " " << event[iPart].daughter2();
151+
LOG(debug) << "Found hadron " << event[iPart].id() << " with rapidity " << rap << " and daughters " << d1 << " " << d2;
124152
// print pdg code of daughters
125153
LOG(debug) << "Daughters: ";
126-
for (int iDau = event[iPart].daughter1(); iDau <= event[iPart].daughter2(); ++iDau) {
154+
for (int iDau = d1; iDau <= d2; ++iDau) {
127155
LOG(debug) << "Daughter " << iDau << ": " << event[iDau].id();
128156
}
129-
bool isCoalDone = CoalescencePythia8(event, mNucleiPdgList, mTrivialCoal, mCoalMomentum, event[iPart].daughter1(), event[iPart].daughter2());
157+
bool isCoalDone = CoalescencePythia8(event, mNucleiPdgList, mTrivialCoal, mCoalMomentum, d1, d2);
130158
if (isCoalDone) {
131-
LOG(debug) << "Coalescence process found for hadron " << event[iPart].id() << " with daughters " << event[iPart].daughter1() << " " << event[iPart].daughter2();
159+
LOG(debug) << "Coalescence process found for hadron " << event[iPart].id() << " with daughters " << d1 << " " << d2;
132160
LOG(debug) << "Check updated daughters: ";
133-
for (int iDau = event[iPart].daughter1(); iDau <= event[iPart].daughter2(); ++iDau) {
161+
for (int iDau = d1; iDau <= d2; ++iDau) {
134162
LOG(debug) << "Daughter " << iDau << ": " << event[iDau].id();
135163
}
136164
return true;
@@ -171,4 +199,4 @@ FairGenerator *generateHFHadToNuclei(int input_trigger_ratio = 5, std::vector<in
171199
myGen->readString("Random:setSeed on");
172200
myGen->readString("Random:seed " + std::to_string(seed));
173201
return myGen;
174-
}
202+
}

0 commit comments

Comments
 (0)