ATLAS Offline Software
Loading...
Searching...
No Matches
ParticleJetTools::FatVertex::FatVertex Struct Reference

#include <FatVertex.h>

Collaboration diagram for ParticleJetTools::FatVertex::FatVertex:

Public Member Functions

 FatVertex (std::vector< const xAOD::TruthParticle * > in_parts, std::vector< const xAOD::TruthParticle * > intern, std::vector< const xAOD::TruthParticle * > out_parts, bool pv)
float getPT () const
void doParentChildLinks (std::vector< FatVertex > &vertices)
int getNumChargedDecays (float minPt=0.0) const
int getNumNeutralDecays (float minPt=0.0, bool include_neutrinos=false) const
VertexType getType (const SG::ConstAccessor< int > &uid_accessor) const
bool operator< (const FatVertex &other) const
bool operator== (const FatVertex &other) const

Public Attributes

std::vector< const xAOD::TruthParticle * > inparts
std::vector< const xAOD::TruthParticle * > internal
std::vector< const xAOD::TruthParticle * > outparts
bool is_pv
const FatVertexparent =nullptr
std::vector< const FatVertex * > children

Private Member Functions

std::optional< DetailedVertexTypeclassifyTauVertex (DetailedVertexType parent_type, bool has_parent) const
std::optional< DetailedVertexTypeclassifyStrangeVertex (DetailedVertexType parent_type, bool has_parent) const
std::optional< DetailedVertexTypeclassifyPhotonVertex (size_t truth_out_parts, size_t truth_valid_out_parts, const std::vector< const xAOD::TruthParticle * > &truth_children) const
std::optional< DetailedVertexTypeclassifyLeptonVertex (const SG::ConstAccessor< int > &uid_accessor, size_t truth_out_parts, size_t truth_valid_out_parts, const std::vector< const xAOD::TruthParticle * > &truth_children) const

Detailed Description

Definition at line 120 of file FatVertex.h.

Constructor & Destructor Documentation

◆ FatVertex()

ParticleJetTools::FatVertex::FatVertex::FatVertex ( std::vector< const xAOD::TruthParticle * > in_parts,
std::vector< const xAOD::TruthParticle * > intern,
std::vector< const xAOD::TruthParticle * > out_parts,
bool pv )
inline

Definition at line 128 of file FatVertex.h.

128 :
129 inparts(std::move(in_parts)), internal(std::move(intern)),
130 outparts(std::move(out_parts)), is_pv(pv) {}
std::vector< const xAOD::TruthParticle * > inparts
Definition FatVertex.h:121
std::vector< const xAOD::TruthParticle * > internal
Definition FatVertex.h:122
std::vector< const xAOD::TruthParticle * > outparts
Definition FatVertex.h:123

Member Function Documentation

◆ classifyLeptonVertex()

std::optional< DetailedVertexType > ParticleJetTools::FatVertex::FatVertex::classifyLeptonVertex ( const SG::ConstAccessor< int > & uid_accessor,
size_t truth_out_parts,
size_t truth_valid_out_parts,
const std::vector< const xAOD::TruthParticle * > & truth_children ) const
private

Definition at line 495 of file FatVertex.cxx.

499 {
500 // Standard lep->lep+gamma : Bremsstrahlung
501 if ( outparts.size() == 2 &&
502 ((outparts[0]->isPhoton() && outparts[1]->isChLepton()) ||
503 (outparts[1]->isPhoton() && outparts[0]->isChLepton())) ) {
505 }
506 // e->e+gamma but the electron doesn't stay in the truth record
507 if ( (truth_out_parts == 2 || truth_out_parts == 1) && truth_valid_out_parts == 1 &&
508 outparts[0]->isPhoton() ) {
510 }
511 // l -> gamma + l which is still brem
512 if ( truth_valid_out_parts == 2 &&
513 ((truth_children[0]->isPhoton() && truth_children[1]->pdgId() == inparts[0]->pdgId()) ||
514 (truth_children[1]->isPhoton() && truth_children[0]->pdgId() == inparts[0]->pdgId())) ) {
516 }
517 // e+ e- annihilation: e+ -> gamma gamma
518 if ( inparts[0]->pdgId() == -11 && outparts.size() == 2 &&
519 outparts[0]->isPhoton() && outparts[1]->isPhoton() ) {
521 }
522 // e+ -> gamma (boosted, partner low energy)
523 if ( inparts[0]->pdgId() == -11 && outparts.size() == 1 &&
524 outparts[0]->isPhoton() ) {
526 }
527 // e+ -> child is single photon
528 if ( inparts[0]->pdgId() == -11 && inparts[0]->nChildren() == 1 &&
529 inparts[0]->child(0)->isPhoton() ) {
531 }
532 // Lepton scattering via Geant (same base UID → material interaction)
533 int in_base_bar = uid_accessor(*inparts[0]) % HepMC::SIM_REGENERATION_INCREMENT;
534 if ( truth_valid_out_parts == 2 && truth_out_parts == 2 &&
535 (in_base_bar == uid_accessor(*outparts[0]) % HepMC::SIM_REGENERATION_INCREMENT ||
536 in_base_bar == uid_accessor(*outparts[1]) % HepMC::SIM_REGENERATION_INCREMENT) ) {
538 }
539 return std::nullopt;
540 }
bool isPhoton(const T &p)
Definition AtlasPID.h:376
constexpr int SIM_REGENERATION_INCREMENT
Constant defining the barcode threshold for regenerated particles, i.e. particles surviving an intera...

◆ classifyPhotonVertex()

std::optional< DetailedVertexType > ParticleJetTools::FatVertex::FatVertex::classifyPhotonVertex ( size_t truth_out_parts,
size_t truth_valid_out_parts,
const std::vector< const xAOD::TruthParticle * > & truth_children ) const
private

Definition at line 452 of file FatVertex.cxx.

455 {
456 // Gamma -> e+e- : Conversion
457 if ( outparts.size() == 2 && outparts[0]->isChLepton() && outparts[1]->isChLepton() &&
458 (outparts[0]->pdgId() + outparts[1]->pdgId()) == 0 ) {
460 }
461 // Conversion from truth children (handles cases where fat vertex merging changed outparts)
462 if ( truth_valid_out_parts == 2 &&
463 (truth_children[0]->pdgId() + truth_children[1]->pdgId()) == 0 ) {
465 }
466 // Photon -> single electron: photoelectric emission
467 if ( std::abs(outparts[0]->pdgId()) == 11 ) {
469 }
470 // Conversion where one lepton was absorbed by material
471 if ( outparts.size() == 1 && internal.size() == 1 &&
472 internal[0]->isChLepton() && outparts[0]->isChLepton() ) {
474 }
475 // Conversion with missing truth record for one lepton
476 if ( outparts.size() == 1 && internal.empty() && truth_out_parts == 2 &&
477 truth_valid_out_parts == 1 && outparts[0]->isChLepton() ) {
479 }
480 // gamma -> gamma + e : Compton scattering
481 if ( outparts.size() == 2 &&
482 ((outparts[0]->isPhoton() && outparts[1]->pdgId() == 11) ||
483 (outparts[1]->isPhoton() && outparts[0]->pdgId() == 11)) ) {
485 }
486 // Catch-all for photon vertices with 3+ outparts containing electrons:
487 // conversion + bremsstrahlung merged by fat vertex builder
488 for (auto p : outparts) {
489 if (p && p->isChLepton()) return DetailedVertexType::Conversion;
490 }
491 // Photon interaction with no electrons — generic material interaction
493 }

◆ classifyStrangeVertex()

std::optional< DetailedVertexType > ParticleJetTools::FatVertex::FatVertex::classifyStrangeVertex ( DetailedVertexType parent_type,
bool has_parent ) const
private

Definition at line 430 of file FatVertex.cxx.

432 {
433 if ( outparts.size() == 1 && outparts[0]->isStrangeHadron() )
435 if ( !has_parent ) return DetailedVertexType::StrangeDecay;
436 switch (parent_type) {
447 default:
449 }
450 }

◆ classifyTauVertex()

std::optional< DetailedVertexType > ParticleJetTools::FatVertex::FatVertex::classifyTauVertex ( DetailedVertexType parent_type,
bool has_parent ) const
private

Definition at line 411 of file FatVertex.cxx.

413 {
414 int num_child_taus = 0;
415 int num_child_photons = 0;
416 for ( auto child : outparts ) {
417 if ( child->isTau() ) num_child_taus++;
418 else if ( child->isPhoton() ) num_child_photons++;
419 }
420 if ( num_child_taus == 0 ) {
421 if (has_parent && parent_type == DetailedVertexType::BHadronDecay)
424 }
425 if ( num_child_taus == 1 && num_child_photons == 1 )
427 return std::nullopt;
428 }

◆ doParentChildLinks()

void ParticleJetTools::FatVertex::FatVertex::doParentChildLinks ( std::vector< FatVertex > & vertices)

Definition at line 558 of file FatVertex.cxx.

558 {
559 // Iterates through all the vertices to set parent and children
560
561 FatVertex* found_parent = nullptr;
562 std::vector<const FatVertex*> found_children;
563
564 for (auto& vertex2 : vertices) {
565 if ( this == &vertex2) continue; // Compare addresses, not objects
566
567 // Check if any of vertex2's outputs match vertex's inputs
568 if (std::find(vertex2.outparts.begin(), vertex2.outparts.end(), inparts[0]) != vertex2.outparts.end()) {
569 if (found_parent) {
570 throw std::runtime_error("Multiple parents detected for a vertex");
571 }
572 found_parent = &vertex2;
573 }
574
575 // Check if vertex2 is a child
576 if (std::find(outparts.begin(), outparts.end(), vertex2.inparts[0]) != outparts.end()) {
577 found_children.push_back(&vertex2);
578 }
579 }
580
581 parent = found_parent;
582 children = std::move(found_children);
583 }
FatVertex(std::vector< const xAOD::TruthParticle * > in_parts, std::vector< const xAOD::TruthParticle * > intern, std::vector< const xAOD::TruthParticle * > out_parts, bool pv)
Definition FatVertex.h:128
std::vector< const FatVertex * > children
Definition FatVertex.h:126

◆ getNumChargedDecays()

int ParticleJetTools::FatVertex::FatVertex::getNumChargedDecays ( float minPt = 0.0) const

Definition at line 285 of file FatVertex.cxx.

285 {
286 int num_charged = 0;
287 for ( auto part : outparts ) {
288 if ( part->isCharged() && part->pt() >= minPt ) num_charged++;
289 }
290 return num_charged;
291 }

◆ getNumNeutralDecays()

int ParticleJetTools::FatVertex::FatVertex::getNumNeutralDecays ( float minPt = 0.0,
bool include_neutrinos = false ) const

Definition at line 293 of file FatVertex.cxx.

293 {
294 int num_neutral = 0;
295 for ( auto part : outparts ) {
296 if ( part->isNeutral() && part->pt() >= minPt && (include_neutrinos || !part->isNeutrino()) )
297 num_neutral++;
298 }
299 return num_neutral;
300 }

◆ getPT()

float ParticleJetTools::FatVertex::FatVertex::getPT ( ) const

Definition at line 280 of file FatVertex.cxx.

280 {
281 if ( inparts.size() == 1 ) { return inparts[0]->pt(); }
282 return sum_4vec(inparts).Pt();
283 }
TLorentzVector sum_4vec(const std::vector< const xAOD::TruthParticle * > &parts)
Definition FatVertex.cxx:14

◆ getType()

VertexType ParticleJetTools::FatVertex::FatVertex::getType ( const SG::ConstAccessor< int > & uid_accessor) const

Definition at line 302 of file FatVertex.cxx.

304 {
306 if ( inparts.size() != 1 ) return DetailedVertexType::OtherNInparts;
307 if ( outparts.empty() ) return DetailedVertexType::OtherNoOutparts;
308
309 // Get first valid child
310 const xAOD::TruthParticle* input_child = nullptr;
311 for (auto child : inparts[0]->decayVtx()->particles_out()) {
312 if (child) { input_child = child; break; }
313 }
314 if ( !input_child ) return DetailedVertexType::Other;
315
316 bool has_parent = parent != nullptr;
317 DetailedVertexType parent_type = has_parent
318 ? parent->getType(uid_accessor).detailedType
320
321 // Heavy flavor: B and C hadron decays
323 if ( inparts[0]->isCharmHadron() ) {
324 if (has_parent && parent_type == DetailedVertexType::BHadronDecay)
327 }
328
329 // Tau decay (can fall through if it's tau bremsstrahlung that doesn't match)
330 if ( inparts[0]->isTau() ) {
331 if (auto t = classifyTauVertex(parent_type, has_parent)) return *t;
332 }
333
334 // Compute shared context for remaining classifiers
335 size_t truth_out_parts = inparts[0]->nChildren();
336 size_t truth_valid_out_parts = num_valid_children(inparts[0]);
337 std::vector<const xAOD::TruthParticle*> truth_children;
338 for ( auto child : inparts[0]->decayVtx()->particles_out() ) {
339 if ( child ) truth_children.push_back(child);
340 }
341 TLorentzVector in4 = inparts[0]->p4();
342 TLorentzVector out4 = sum_4vec(inparts[0]->decayVtx()->particles_out());
343 float deltaE = std::abs(in4.E() - out4.E());
344
345 // Energy-based material interaction (all valid children, large energy loss)
346 if ( truth_out_parts == truth_valid_out_parts && truth_valid_out_parts > 1 ) {
347 if ( deltaE > 100.0 ) return DetailedVertexType::MaterialInteraction;
348 }
349
350 // Mass-based material interaction: if the sum of children rest masses
351 // exceeds the parent rest mass, mass was created from kinetic energy
352 // via nuclear interaction with detector material. In a genuine decay,
353 // parent >= sum(m_children) always holds by conservation.
354 if ( truth_valid_out_parts > 1 ) {
355 double children_mass_sum = 0.0;
356 for (const auto& child : truth_children) {
357 children_mass_sum += child->p4().M();
358 }
359 if ( children_mass_sum > in4.M() + 1.0 ) { // 1 MeV tolerance
361 }
362 }
363
364 // Strange hadron
365 if ( inparts[0]->isStrangeHadron() ) {
366 if (auto t = classifyStrangeVertex(parent_type, has_parent)) return *t;
367 }
368
369 // Photon (conversion, photoelectric, compton)
370 if ( inparts[0]->isPhoton() ) {
371 if (auto t = classifyPhotonVertex(truth_out_parts, truth_valid_out_parts, truth_children)) return *t;
372 }
373
374 // Charged lepton (bremsstrahlung, e+e- annihilation, material interaction)
375 if ( inparts[0]->isChLepton() ) {
376 if (auto t = classifyLeptonVertex(uid_accessor, truth_out_parts, truth_valid_out_parts, truth_children)) return *t;
377 }
378
379 // Pion decay: pi0 -> gamma gamma is always a decay.
380 // Charged pi -> mu + nu: check for muon+neutrino signature.
381 // Anything else from a charged pion (hadronic products, nuclear
382 // fragments) is a material interaction.
383 if ( std::abs(inparts[0]->pdgId()) == 111 ) {
385 }
386 if ( std::abs(inparts[0]->pdgId()) == 211 ) {
387 bool has_muon = false, has_neutrino = false;
388 for (auto c : truth_children) {
389 if (std::abs(c->pdgId()) == 13) has_muon = true;
390 if (std::abs(c->pdgId()) == 14) has_neutrino = true;
391 }
392 if (has_muon && has_neutrino) return DetailedVertexType::PionDecay;
394 }
395
396 // Fallback: simulation particle with missing children
397 if ( truth_out_parts > truth_valid_out_parts && truth_out_parts > 1 &&
398 inparts[0]->child(0) && HepMC::is_simulation_particle(inparts[0]->child(0)) ) {
400 }
401
402 // Single output particle fallbacks
403 if ( truth_out_parts == 1 && truth_valid_out_parts == 1 ) {
404 if ( deltaE > 10 ) return DetailedVertexType::MaybeMaterialInteraction;
406 }
407
409 }
bool isTau(const T &p)
Definition AtlasPID.h:208
bool isChLepton(const T &p)
APID: the fourth generation leptons are leptons.
Definition AtlasPID.h:199
bool isStrangeHadron(const T &p)
Definition AtlasPID.h:910
bool isBottomHadron(const T &p)
Definition AtlasPID.h:912
bool isCharmHadron(const T &p)
Definition AtlasPID.h:911
bool is_simulation_particle(const T &p)
Method to establish if a particle (or barcode) was created during the simulation (TODO update to be s...
int num_valid_children(const xAOD::TruthParticle *part)
Definition FatVertex.cxx:23
TruthParticle_v1 TruthParticle
Typedef to implementation.
std::optional< DetailedVertexType > classifyPhotonVertex(size_t truth_out_parts, size_t truth_valid_out_parts, const std::vector< const xAOD::TruthParticle * > &truth_children) const
std::optional< DetailedVertexType > classifyLeptonVertex(const SG::ConstAccessor< int > &uid_accessor, size_t truth_out_parts, size_t truth_valid_out_parts, const std::vector< const xAOD::TruthParticle * > &truth_children) const
std::optional< DetailedVertexType > classifyStrangeVertex(DetailedVertexType parent_type, bool has_parent) const
std::optional< DetailedVertexType > classifyTauVertex(DetailedVertexType parent_type, bool has_parent) const

◆ operator<()

bool ParticleJetTools::FatVertex::FatVertex::operator< ( const FatVertex & other) const

Definition at line 542 of file FatVertex.cxx.

542 {
543 if ( inparts.size() != other.inparts.size() )
544 return inparts.size() < other.inparts.size();
545 if ( internal.size() != other.internal.size() )
546 return internal.size() < other.internal.size();
547 if ( outparts.size() != other.outparts.size() )
548 return outparts.size() < other.outparts.size();
549 return getPT() < other.getPT();
550 }

◆ operator==()

bool ParticleJetTools::FatVertex::FatVertex::operator== ( const FatVertex & other) const

Definition at line 552 of file FatVertex.cxx.

552 {
553 return this->inparts == other.inparts && this->internal == other.internal &&
554 this->outparts == other.outparts;
555 }

Member Data Documentation

◆ children

std::vector<const FatVertex*> ParticleJetTools::FatVertex::FatVertex::children

Definition at line 126 of file FatVertex.h.

◆ inparts

std::vector<const xAOD::TruthParticle*> ParticleJetTools::FatVertex::FatVertex::inparts

Definition at line 121 of file FatVertex.h.

◆ internal

std::vector<const xAOD::TruthParticle*> ParticleJetTools::FatVertex::FatVertex::internal

Definition at line 122 of file FatVertex.h.

◆ is_pv

bool ParticleJetTools::FatVertex::FatVertex::is_pv

Definition at line 124 of file FatVertex.h.

◆ outparts

std::vector<const xAOD::TruthParticle*> ParticleJetTools::FatVertex::FatVertex::outparts

Definition at line 123 of file FatVertex.h.

◆ parent

const FatVertex* ParticleJetTools::FatVertex::FatVertex::parent =nullptr

Definition at line 125 of file FatVertex.h.


The documentation for this struct was generated from the following files: