ATLAS Offline Software
Reconstruction/tauRecTools/Root/HelperFunctions.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 
7 #include <TObjString.h>
8 #include <TObjArray.h>
9 #include <TFile.h>
10 #include <TTree.h>
11 
12 #include <iostream>
13 
14 namespace tauRecTools {
15  ANA_MSG_SOURCE(msgHelperFunction, "HelperFunction")
16 }
17 
18 
19 
21  const xAOD::Vertex* jetVertex = nullptr;
22 
23  bool isAvailable = jet.getAssociatedObject("OriginVertex", jetVertex);
24  if (! isAvailable) {
25  return nullptr;
26  }
27 
28  return jetVertex;
29 }
30 
31 
32 
33 TLorentzVector tauRecTools::getTauAxis(const xAOD::TauJet& tau, bool doVertexCorrection) {
34  TLorentzVector tauAxis;
35  if (doVertexCorrection) {
37  }
38  else {
40  }
41 
42  return tauAxis;
43 }
44 
45 //________________________________________________________________________________
47  const int flagsize=sizeof(flag)*8;
50  return flag;
51 }
52 
53 //________________________________________________________________________________
55 {
56  //should we be safe and ask if the links are available?
57  const xAOD::TauTrack* xTrack1 = *l1;
58  const xAOD::TauTrack* xTrack2 = *l2;
59 
60  //return classified charged, then isolation (wide tracks), else by pt
63 
64  if(f1==f2)
65  return xTrack1->pt()>xTrack2->pt();
66  return f1<f2;
67 }
68 
69 //________________________________________________________________________________
70 TLorentzVector tauRecTools::GetConstituentP4(const xAOD::JetConstituent& constituent) {
71  using namespace tauRecTools::msgHelperFunction;
72 
73  TLorentzVector constituentP4;
74 
75  if( constituent->type() == xAOD::Type::CaloCluster ) {
76  const xAOD::CaloCluster* cluster = static_cast<const xAOD::CaloCluster*>( constituent->rawConstituent() );
77  constituentP4 = cluster->p4();
78  }
79  else if ( constituent->type() == xAOD::Type::FlowElement ) {
80  const xAOD::FlowElement* fe = static_cast<const xAOD::FlowElement*>( constituent->rawConstituent() );
81  constituentP4 = fe->p4();
82  }
83  else {
84  ANA_MSG_ERROR("GetConstituentP4: Seed jet constituent type not supported!");
85  constituentP4.SetPtEtaPhiE(constituent.pt(), constituent.eta(), constituent.phi(), constituent.e());
86  }
87 
88  return constituentP4;
89 }
90 
91 
92 
94  // build pi0s and shots for 0-5p taus
95  return ( tau.nTracks() <= 5 );
96 }
tauRecTools::getJetVertex
const xAOD::Vertex * getJetVertex(const xAOD::Jet &jet)
Return the vertex of jet candidate.
Definition: Reconstruction/tauRecTools/Root/HelperFunctions.cxx:20
xAOD::TauJetParameters::IntermediateAxis
@ IntermediateAxis
Definition: TauDefs.h:338
xAOD::TauJetParameters::classifiedFake
@ classifiedFake
Definition: TauDefs.h:409
tauRecTools::getTauAxis
TLorentzVector getTauAxis(const xAOD::TauJet &tau, bool doVertexCorrection=true)
Return the four momentum of the tau axis The tau axis is widely used to select clusters and cells in ...
Definition: Reconstruction/tauRecTools/Root/HelperFunctions.cxx:33
tauRecTools::GetConstituentP4
TLorentzVector GetConstituentP4(const xAOD::JetConstituent &constituent)
Definition: Reconstruction/tauRecTools/Root/HelperFunctions.cxx:70
xAOD::TauJet_v3::nTracks
size_t nTracks(TauJetParameters::TauTrackFlag flag=TauJetParameters::TauTrackFlag::classifiedCharged) const
Definition: TauJet_v3.cxx:526
xAOD::TauJetParameters::classifiedCharged
@ classifiedCharged
Definition: TauDefs.h:406
ANA_MSG_ERROR
#define ANA_MSG_ERROR(xmsg)
Macro printing error messages.
Definition: Control/AthToolSupport/AsgMessaging/AsgMessaging/MessageCheck.h:294
xAOD::CaloCluster
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
Definition: Event/xAOD/xAODCaloEvent/xAODCaloEvent/CaloCluster.h:19
ExpressionParsing::isAvailable
bool isAvailable(const T_Aux &cont, SG::auxid_t auxid)
xAOD::JetConstituent::phi
double phi() const
The azimuthal angle ( ) of the particle.
Definition: JetConstituentVector.h:72
read_hist_ntuple.f2
f2
Definition: read_hist_ntuple.py:20
xAOD::CaloCluster_v1
Description of a calorimeter cluster.
Definition: CaloCluster_v1.h:59
skel.l2
l2
Definition: skel.GENtoEVGEN.py:426
jet
Definition: JetCalibTools_PlotJESFactors.cxx:23
xAOD::TauTrack_v1::flagSet
TrackFlagType flagSet() const
xAOD::FlowElement
FlowElement_v1 FlowElement
Definition of the current "pfo version".
Definition: FlowElement.h:16
tauRecTools::doPi0andShots
bool doPi0andShots(const xAOD::TauJet &tau)
Determines whether pi0s and shots should be built for a tau candidate.
Definition: Reconstruction/tauRecTools/Root/HelperFunctions.cxx:93
xAOD::TauJet_v3
Class describing a tau jet.
Definition: TauJet_v3.h:41
master.flag
bool flag
Definition: master.py:29
xAOD::JetConstituent::pt
double pt() const
Definition: JetConstituentVector.h:68
xAOD::JetConstituent::rawConstituent
const IParticle * rawConstituent() const
Access the real underlying IParticle.
Definition: JetConstituentVector.h:102
xAOD::CaloCluster_v1::p4
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition: CaloCluster_v1.cxx:465
ANA_MSG_SOURCE
#define ANA_MSG_SOURCE(NAME, TITLE)
the source code part of ANA_MSG_SOURCE
Definition: Control/AthToolSupport/AsgMessaging/AsgMessaging/MessageCheck.h:133
xAOD::TauTrack_v1::pt
virtual double pt() const
The transverse momentum ( ) of the particle.
xAOD::JetConstituent::type
Type::ObjectType type() const
The full 4-momentum of the particle.
Definition: JetConstituentVector.h:91
xAOD::TauTrack_v1::TrackFlagType
uint16_t TrackFlagType
Definition: TauTrack_v1.h:61
tauRecTools::isolateClassifiedBits
xAOD::TauTrack::TrackFlagType isolateClassifiedBits(xAOD::TauTrack::TrackFlagType flag)
Definition: Reconstruction/tauRecTools/Root/HelperFunctions.cxx:46
xAOD::Jet_v1
Class describing a jet.
Definition: Jet_v1.h:57
xAOD::FlowElement_v1::p4
virtual FourMom_t p4() const override
The full 4-momentum of the particle.
Definition: FlowElement_v1.cxx:33
xAOD::Vertex_v1
Class describing a Vertex.
Definition: Vertex_v1.h:42
xAOD::TauTrack_v1
Definition: TauTrack_v1.h:27
xAOD::TauJet_v3::p4
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition: TauJet_v3.cxx:97
HelperFunctions.h
xAOD::JetConstituent::e
double e() const
The total energy of the particle.
Definition: JetConstituentVector.h:76
xAOD::JetConstituent
4-vector of jet constituent at the scale used during jet finding.
Definition: JetConstituentVector.h:61
xAOD::JetConstituent::eta
double eta() const
The pseudorapidity ( ) of the particle.
Definition: JetConstituentVector.h:70
skel.l1
l1
Definition: skel.GENtoEVGEN.py:425
tauRecTools
Implementation of a TrackClassifier based on an RNN.
Definition: BDTHelper.cxx:12
xAOD::TauJetParameters::DetectorAxis
@ DetectorAxis
Definition: TauDefs.h:337
read_hist_ntuple.f1
f1
Definition: read_hist_ntuple.py:4
tauRecTools::sortTracks
bool sortTracks(const ElementLink< xAOD::TauTrackContainer > &l1, const ElementLink< xAOD::TauTrackContainer > &l2)
Definition: Reconstruction/tauRecTools/Root/HelperFunctions.cxx:54
xAOD::FlowElement_v1
A detector object made of other lower level object(s)
Definition: FlowElement_v1.h:25