ATLAS Offline Software
Loading...
Searching...
No Matches
ParticleJetTools Namespace Reference

Classes

class  IParticleLinker
struct  LabelDecorators
struct  LabelNames
struct  Particles

Functions

Amg::Vector3D p3 (const xAOD::TruthVertex *p)
Amg::Vector3D p3 (const xAOD::Jet &j)
Amg::Vector3D signalProcessP3 (const xAOD::TruthEventContainer &)
void setJetLabels (const xAOD::Jet &jet, const Particles &particles, const LabelNames &names)
void setJetLabels (const xAOD::Jet &jet, const Particles &particles, const LabelDecorators &decs)
float partPt (const xAOD::TruthParticle *part)
float partLxy (const xAOD::TruthParticle *part, const Amg::Vector3D &origin)
float partDR (const xAOD::TruthParticle *part, const xAOD::Jet &jet)
int partPdgId (const xAOD::TruthParticle *part)
float positionDPhi (const xAOD::TruthParticle *part, const xAOD::Jet &jet, const Amg::Vector3D &origin)
float positionDEta (const xAOD::TruthParticle *part, const xAOD::Jet &jet, const Amg::Vector3D &origin)
void childrenRemoved (const std::vector< const xAOD::TruthParticle * > &parents, std::vector< const xAOD::TruthParticle * > &children)
template<typename T>
void declareProperties (T &tool, LabelNames *n)
bool isChild (const xAOD::TruthParticle *p, const xAOD::TruthParticle *c)

Function Documentation

◆ childrenRemoved()

void ParticleJetTools::childrenRemoved ( const std::vector< const xAOD::TruthParticle * > & parents,
std::vector< const xAOD::TruthParticle * > & children )

Definition at line 63 of file ParticleJetLabelCommon.cxx.

66 {
67
68 if ( &parents == &children ) {
69 // Same vector provided for both inputs. Extra care needed in
70 // this case...
71 const std::vector<const xAOD::TruthParticle*> copyParents = parents;
72 childrenRemoved(copyParents, children);
73 return;
74 }
75 // We can now safely assume that parents will not be modified
76 // during this loop.
77 for ( const xAOD::TruthParticle* p : parents ) {
78 if (!p) continue;
79 children.erase(std::remove_if(children.begin(),
80 children.end(),
81 [p](const xAOD::TruthParticle* c) { return (c && isChild(p, c)); }),
82 children.end());
83 // auto erased = std::erase_if(children, [p](const xAOD::TruthParticle* c) { return (c && isChild(p, c)); }); // C++20
84 }
85 return;
86 }
void childrenRemoved(const std::vector< const xAOD::TruthParticle * > &parents, std::vector< const xAOD::TruthParticle * > &children)
DataModel_detail::iterator< DVL > remove_if(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end, Predicate pred)
Specialization of remove_if for DataVector/List.
TruthParticle_v1 TruthParticle
Typedef to implementation.

◆ declareProperties()

template<typename T>
void ParticleJetTools::declareProperties ( T & tool,
LabelNames * n )

Definition at line 109 of file ParticleJetLabelCommon.h.

109 {
110 tool.declareProperty("LabelName", n->singleint="", "Jet label attribute to be added.");
111 tool.declareProperty("DoubleLabelName", n->doubleint="", "Jet label attribute to be added (with the possibility of up to 2 matched hadrons).");
112 tool.declareProperty("LabelPtName", n->pt="", "Attribute for labelling particle pt");
113 tool.declareProperty("LabelLxyName", n->Lxy="", "Attribute for Lxy of labelling particle");
114 tool.declareProperty("LabelDRName", n->dr="", "Attribute for dR(part, jet) for labelling particle");
115 tool.declareProperty("LabelPdgIdName", n->pdgId="", "Attribute for pdgID of labelling particle");
116 tool.declareProperty("LabelPositionDPhiName", n->positionDPhi="", "Attribute for the position dPhi of the labeling particle ");
117 tool.declareProperty("LabelPositionDEtaName", n->positionDEta="", "Attribute for the position dEta of the labeling particle ");
118 tool.declareProperty("LabelBarcodeName", n->uniqueID="", "Attribute for uniqueID of labeling particle");
119 tool.declareProperty("ChildLxyName", n->childLxy="", "Attribute for the labeling particle child Lxy");
120 tool.declareProperty("ChildPtName", n->childPt="", "Attribute for the labeling particle child Pt");
121 tool.declareProperty("ChildPdgIdName", n->childPdgId="", "Attribute for the labeling particle child pdg ID");
122 tool.declareProperty("ChildPositionDPhiName", n->childPositionDPhi="", "Attribute for the position dPhi of the labeling particle child");
123 tool.declareProperty("ChildPositionDEtaName", n->childPositionDEta="", "Attribute for the position dEta of the labeling particle child");
124 // ATLASRECTS-8290: this is for backward compatability, remove eventually
125 tool.declareProperty("useBarcode", n->useBarcode=false, "use barcode instead of uid");
126 }

◆ isChild()

bool ParticleJetTools::isChild ( const xAOD::TruthParticle * p,
const xAOD::TruthParticle * c )
inline

Definition at line 34 of file ParticleJetLabelCommon.cxx.

34 {
35
36 if ( HepMC::is_same_particle(p,c) ) { return false; }
37
38 for (size_t iC = 0; iC < p->nChildren(); iC++) {
39 const xAOD::TruthParticle* cc = p->child(iC);
40 if (!cc) { continue; }
41
42 if ( HepMC::is_same_particle(cc,c) ) { return true; }
43
44 if (isChild(cc, c)) { return true; }
45 }
46
47 return false;
48 }
bool is_same_particle(const T1 &p1, const T2 &p2)
Method to establish if two particles in the GenEvent actually represent the same particle.
bool isChild(const xAOD::TruthParticle *p, const xAOD::TruthParticle *c)

◆ p3() [1/2]

Amg::Vector3D ParticleJetTools::p3 ( const xAOD::Jet & j)

Definition at line 54 of file ParticleJetLabelCommon.cxx.

54 {
55 return {j.px(), j.py(), j.pz()};
56 }
float py() const
The y-component of the jet's momentum.
Definition Jet_v1.cxx:94
float px() const
The x-component of the jet's momentum.
Definition Jet_v1.cxx:90
float pz() const
The z-component of the jet's momentum.
Definition Jet_v1.cxx:99

◆ p3() [2/2]

Amg::Vector3D ParticleJetTools::p3 ( const xAOD::TruthVertex * p)

Definition at line 50 of file ParticleJetLabelCommon.cxx.

50 {
51 if (!p) return {NAN, NAN, NAN};
52 return {p->x(), p->y(), p->z()};
53 }

◆ partDR()

float ParticleJetTools::partDR ( const xAOD::TruthParticle * part,
const xAOD::Jet & jet )

Definition at line 298 of file ParticleJetLabelCommon.cxx.

298 {
299 if (!part) return NAN;
300 return part->p4().DeltaR(jet.p4());
301 }

◆ partLxy()

float ParticleJetTools::partLxy ( const xAOD::TruthParticle * part,
const Amg::Vector3D & origin )

Definition at line 293 of file ParticleJetLabelCommon.cxx.

293 {
294 if (!part) return NAN;
295 if (const auto* vx = part->decayVtx() ) return (p3(vx) - origin).perp();
296 return INFINITY;
297 }
Amg::Vector3D p3(const xAOD::TruthVertex *p)

◆ partPdgId()

int ParticleJetTools::partPdgId ( const xAOD::TruthParticle * part)

Definition at line 302 of file ParticleJetLabelCommon.cxx.

302 {
303 if (!part) return 0;
304 return part->pdgId();
305 }

◆ partPt()

float ParticleJetTools::partPt ( const xAOD::TruthParticle * part)

Definition at line 289 of file ParticleJetLabelCommon.cxx.

289 {
290 if (!part) return NAN;
291 return part->pt();
292 }

◆ positionDEta()

float ParticleJetTools::positionDEta ( const xAOD::TruthParticle * part,
const xAOD::Jet & jet,
const Amg::Vector3D & origin )

Definition at line 316 of file ParticleJetLabelCommon.cxx.

318 {
319 if (!part) return NAN;
320 if (const auto* vx = part->decayVtx() ) {
321 Amg::Vector3D displacement = p3(vx) - origin;
322 return displacement.eta() - jet.eta();
323 }
324 return INFINITY;
325 }
Eigen::Matrix< double, 3, 1 > Vector3D

◆ positionDPhi()

float ParticleJetTools::positionDPhi ( const xAOD::TruthParticle * part,
const xAOD::Jet & jet,
const Amg::Vector3D & origin )

Definition at line 306 of file ParticleJetLabelCommon.cxx.

308 {
309 if (!part) return NAN;
310 if (const auto* vx = part->decayVtx() ) {
311 Amg::Vector3D displacement = p3(vx) - origin;
312 return p3(jet).deltaPhi(displacement);
313 }
314 return INFINITY;
315 }

◆ setJetLabels() [1/2]

void ParticleJetTools::setJetLabels ( const xAOD::Jet & jet,
const Particles & particles,
const LabelDecorators & decs )

Definition at line 164 of file ParticleJetLabelCommon.cxx.

166 {
167
168 // we also want to save information about the maximum pt particle of the labeling partons
169 auto getMaxPtPart = [](const auto& container) -> const xAOD::TruthParticle* {
170 if (container.size() == 0) return nullptr;
171 auto itr = std::max_element(container.begin(), container.end(),
172 [](auto* p1, auto* p2) {
173 return p1->pt() < p2->pt();
174 });
175 return *itr;
176 };
177
178 // set truth label for jets above pt threshold
179 // hierarchy: b > c > tau > light
180 int label = 0; // default: light
181 const xAOD::TruthParticle* labelling_particle = nullptr;
182 const xAOD::TruthParticle* child_particle = nullptr;
183 if (particles.b.size()) {
184 label = 5;
185 labelling_particle = getMaxPtPart(particles.b);
186 child_particle = getCharmChild(labelling_particle);
187 } else if (particles.c.size()) {
188 label = 4;
189 labelling_particle = getMaxPtPart(particles.c);
190 } else if (particles.tau.size()) {
191 label = 15;
192 labelling_particle = getMaxPtPart(particles.tau);
193 }
194
195 const Amg::Vector3D& origin = particles.origin;
196
197 // decorate info about the labelling particle
198 decs.singleint(jet) = label;
199 if (label == 0) {
200 decs.pt(jet) = NAN;
201 decs.Lxy(jet) = NAN;
202 decs.dr(jet) = NAN;
203 decs.pdgId(jet) = 0;
204 decs.positionDPhi(jet) = NAN;
205 decs.positionDEta(jet) = NAN;
207 decs.childLxy(jet) = NAN;
208 decs.childPt(jet) = NAN;
209 decs.childPdgId(jet) = 0;
210 decs.childPositionDPhi(jet) = NAN;
211 decs.childPositionDEta(jet) = NAN;
212 } else {
213 decs.pt(jet) = partPt(labelling_particle);
214 decs.Lxy(jet) = partLxy(labelling_particle, origin);
215 decs.dr(jet) = partDR(labelling_particle, jet);
216 decs.pdgId(jet) = partPdgId(labelling_particle);
217 decs.positionDPhi(jet) = positionDPhi(labelling_particle, jet, origin);
218 decs.positionDEta(jet) = positionDEta(labelling_particle, jet, origin);
219 decs.uniqueID(jet) = labelling_particle ?
220 // ATLASRECTS-8290: this should be replaced with ->uid()
221 decs.acc_uid(*labelling_particle) : HepMC::INVALID_PARTICLE_ID;
222 decs.childLxy(jet) = partLxy(child_particle, origin);
223 decs.childPt(jet) = partPt(child_particle);
224 decs.childPdgId(jet) = partPdgId(child_particle);
225 decs.childPositionDPhi(jet) = positionDPhi(child_particle, jet, origin);
226 decs.childPositionDEta(jet) = positionDEta(child_particle, jet, origin);
227 }
228
229 // extended flavour label
230 int double_label = 0;
231 if (particles.b.size()) {
232 if (particles.b.size() >= 2)
233 double_label = 55;
234
235 else if (particles.c.size())
236 double_label = 54;
237
238 else
239 double_label = 5;
240
241 } else if (particles.c.size()) {
242 if (particles.c.size() >= 2)
243 double_label = 44;
244
245 else
246 double_label = 4;
247
248 } else if (particles.tau.size()){
249
250 bool hasElectrondecay = false;
251 bool hasMuondecay = false;
252 bool hasHadronicdecay = false;
253 for (auto itau: particles.tau){
254 for (size_t i = 0; i< itau->nChildren(); i++){ // tau children loop
255 if (itau->child(i)->absPdgId() == 12 || itau->child(i)->absPdgId() == 14 || itau->child(i)->absPdgId() == 16) continue;
256 if (MC::isMuon(itau->child(i))) hasMuondecay = true;
257 else if (MC::isElectron(itau->child(i))) hasElectrondecay = true;
258 else hasHadronicdecay = true;
259 }
260 }
261
262 if (particles.tau.size() >= 2){
263 // check if we have at least one hadronic tau
264 if (hasHadronicdecay){
265 // check if we have also tau->mu decay
266 if (hasMuondecay) double_label = 151513;
267 // check if we have also tau->el decay
268 else if (hasElectrondecay) double_label = 151511;
269 // otherwise fully hadronic di-tau decay
270 else double_label = 1515;
271 }
272 } else {
273 // only consider hadronic tau decay
274 if (hasHadronicdecay) double_label = 15;
275 }
276
277 }
278
279 decs.doubleint(jet) = double_label;
280
281 }
std::string label(const std::string &format, int i)
Definition label.h:19
constexpr int INVALID_PARTICLE_ID
bool isElectron(const T &p)
bool isMuon(const T &p)
float positionDEta(const xAOD::TruthParticle *part, const xAOD::Jet &jet, const Amg::Vector3D &origin)
float partDR(const xAOD::TruthParticle *part, const xAOD::Jet &jet)
float partPt(const xAOD::TruthParticle *part)
float partLxy(const xAOD::TruthParticle *part, const Amg::Vector3D &origin)
float positionDPhi(const xAOD::TruthParticle *part, const xAOD::Jet &jet, const Amg::Vector3D &origin)
int partPdgId(const xAOD::TruthParticle *part)
SG::AuxElement::Decorator< float > childPositionDEta
SG::AuxElement::Decorator< int > uniqueID
SG::AuxElement::Decorator< int > pdgId
SG::AuxElement::Decorator< float > pt
SG::AuxElement::Decorator< int > doubleint
SG::AuxElement::Decorator< float > Lxy
SG::AuxElement::Decorator< float > dr
SG::AuxElement::Decorator< int > childPdgId
SG::AuxElement::Decorator< float > positionDEta
SG::AuxElement::Decorator< float > childPositionDPhi
SG::AuxElement::Decorator< float > positionDPhi
SG::AuxElement::Decorator< float > childLxy
SG::AuxElement::Decorator< float > childPt
SG::AuxElement::Decorator< int > singleint

◆ setJetLabels() [2/2]

void ParticleJetTools::setJetLabels ( const xAOD::Jet & jet,
const Particles & particles,
const LabelNames & names )

Definition at line 283 of file ParticleJetLabelCommon.cxx.

285 {
286 setJetLabels(jet, particles, LabelDecorators(names));
287 }
void setJetLabels(const xAOD::Jet &jet, const Particles &particles, const LabelNames &names)

◆ signalProcessP3()

Amg::Vector3D ParticleJetTools::signalProcessP3 ( const xAOD::TruthEventContainer & events)

Definition at line 57 of file ParticleJetLabelCommon.cxx.

57 {
58 if (events.empty()) return {NAN, NAN, NAN};
59 return p3(events.at(0)->signalProcessVertex());
60 }