ATLAS Offline Software
Loading...
Searching...
No Matches
GeneralTauPlots.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5
6#include "GeneralTauPlots.h"
8
9namespace Tau{
10
11GeneralTauPlots::GeneralTauPlots(PlotBase* pParent, const std::string& sDir, const std::string& sTauJetContainerName):
12 PlotBase(pParent, sDir),
13 m_oTauKinematicPlots(this, "", sTauJetContainerName),
14 m_sTauJetContainerName(sTauJetContainerName)
15{
16}
17
21
23
24 m_tauCharge = Book1D("Charge",m_sTauJetContainerName + " Tau charge; charge; # Taus",7,-3.,4.);
25 m_tauNChargedTracks = Book1D("NChargedTracks", m_sTauJetContainerName + " Tau n charged tracks; nCharged; # Taus",10,0.,10.);
26 m_tauNIsolatedTracks = Book1D("NIsolatedTracks", m_sTauJetContainerName + " Tau n isolated tracks; nIsolated; # Taus",10,0.,10.);
27 m_tauNCoreTracks = Book1D("NCoreTracks",m_sTauJetContainerName + " Tau n core tracks; nCore; # Taus",15,0.,15.);
28 m_tauNWideTracks = Book1D("NWideTracks",m_sTauJetContainerName + " Tau n wide tracks; nWide; # Taus",20,0.,20.);
29 m_ptHighPt = Book1D("ptHighPt", m_sTauJetContainerName+" HighPt"+"; pt; # Taus",20, 0.0, 1500.0);
30 m_RNNEleScore = Book1D("RNNEleScore", m_sTauJetContainerName+" RNNEleScore;RNNEleScore;# Tau", 50,0.,1.);
31 m_RNNEleScoreSigTrans = Book1D("RNNEleScoreSigTrans", m_sTauJetContainerName+" RNNEleScoreSigTrans;RNNEleScoreSigTrans;"+"# Tau", 50,0.,1.);
32 m_GNTauScore = Book1D("GNTauScore", m_sTauJetContainerName+" GNTauScore;GNTauScore;# Tau", 50,0.,1.);
33 m_GNTauScoreSigTrans = Book1D("GNTauScoreSigTrans", m_sTauJetContainerName+" GNTauScoreSigTrans;GNTauScoreSigTrans;"+"# Tau", 50,0.,1.);
34 m_ptGNTauLoose = Book1D("ptGNTauSigLoose",m_sTauJetContainerName+" GNTauSigLoose; pt; # Taus", 20, 0.0, 150.0);
35 m_ptGNTauLooseHighPt = Book1D("ptGNTauSigLooseHighPt", m_sTauJetContainerName+" GNTauSigLooseHighPt; pt"+"; # Taus",20, 0.0, 1500.0);
36 m_ptGNTauMedium = Book1D("ptGNTauSigMedium",m_sTauJetContainerName+" GNTauSigMedium; pt; # Taus", 20, 0.0, 150.0);
37 m_ptGNTauMediumHighPt = Book1D("ptGNTauSigMediumHighPt", m_sTauJetContainerName+" GNTauSigMediumHighPt; pt"+"; # Taus",20, 0.0, 1500.0);
38 m_ptGNTauTight = Book1D("ptGNTauSigTight",m_sTauJetContainerName+" GNTauSigTight; pt; # Taus", 20, 0.0, 150.0);
39 m_ptGNTauTightHighPt = Book1D("ptGNTauSigTightHighPt", m_sTauJetContainerName+" GNTauSigTightHighPt; pt"+"; # Taus",20, 0.0, 1500.0);
40
41
42
43}
44
45void GeneralTauPlots::fill(const xAOD::TauJet& tau, float weight) {
46 m_oTauKinematicPlots.fill(tau, weight);
47 m_tauCharge->Fill(tau.charge(), weight);
48 m_tauNChargedTracks->Fill(tau.nTracks(), weight);
52 m_ptHighPt->Fill(tau.pt()/Athena::Units::GeV, weight);
53
54 static const SG::ConstAccessor<float> acc_RNNEleScore("RNNEleScore");
55 if ( acc_RNNEleScore.isAvailable(tau) ) {
57 if ( rnnScore > -2.0 ) m_RNNEleScore->Fill(rnnScore, weight);
58 }
59 static const SG::ConstAccessor<float> acc_RNNEleScoreSigTrans("RNNEleScoreSigTrans_v1");
60 if ( acc_RNNEleScoreSigTrans.isAvailable(tau) ) {
61 float rnnScore = acc_RNNEleScoreSigTrans(tau);
62 m_RNNEleScoreSigTrans->Fill(rnnScore, weight);
63 }
64 static const SG::ConstAccessor<float> acc_GNTauScore("GNTauScore_v0prune");
65 if ( acc_GNTauScore.isAvailable(tau) ) {
66 float gntauScore = acc_GNTauScore(tau);
67 if ( gntauScore > -2.0 ) m_GNTauScore->Fill(gntauScore, weight);
68 }
69 static const SG::ConstAccessor<float> acc_GNTauScoreSigTrans("GNTauScoreSigTrans_v0prune");
70 if ( acc_GNTauScoreSigTrans.isAvailable(tau) ) {
71 float gntauScoreSigTrans = acc_GNTauScoreSigTrans(tau);
72 if ( gntauScoreSigTrans > -2.0 ) m_GNTauScoreSigTrans->Fill(gntauScoreSigTrans, weight);
73 }
74
75 static const SG::ConstAccessor<char> acc_GNTauL("GNTauL_v0prune");
76 if( acc_GNTauL.isAvailable(tau) && acc_GNTauL(tau)) {
77 m_ptGNTauLoose ->Fill(tau.pt()/Athena::Units::GeV, weight);
78 m_ptGNTauLooseHighPt->Fill(tau.pt()/Athena::Units::GeV, weight);
79 }
80 static const SG::ConstAccessor<char> acc_GNTauM("GNTauM_v0prune");
81 if( acc_GNTauM.isAvailable(tau) && acc_GNTauM(tau)) {
82 m_ptGNTauMedium ->Fill(tau.pt()/Athena::Units::GeV, weight);
83 m_ptGNTauMediumHighPt->Fill(tau.pt()/Athena::Units::GeV, weight);
84 }
85 static const SG::ConstAccessor<char> acc_GNTauT("GNTauT_v0prune");
86 if( acc_GNTauT.isAvailable(tau) && acc_GNTauT(tau)) {
87 m_ptGNTauTight ->Fill(tau.pt()/Athena::Units::GeV, weight);
88 m_ptGNTauTightHighPt->Fill(tau.pt()/Athena::Units::GeV, weight);
89 }
90}
91
92
93}
Wrapper to avoid constant divisions when using units.
TH1D * Book1D(const std::string &name, const std::string &labels, int nBins, float start, float end, bool prependDir=true)
Book a TH1D histogram.
Definition PlotBase.cxx:94
PlotBase(PlotBase *parent, const std::string &sDir)
Definition PlotBase.cxx:29
Helper class to provide constant type-safe access to aux data.
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
Tau::TauKinematicPlots m_oTauKinematicPlots
std::string m_sTauJetContainerName
void fill(const xAOD::TauJet &tau, float weight)
GeneralTauPlots(PlotBase *pParent, const std::string &sDir, const std::string &sTauJetContainerName)
float charge() const
virtual double pt() const
The transverse momentum ( ) of the particle.
double discriminant(TauJetParameters::TauID discID) const
Get value of discriminant.
size_t nTracks(TauJetParameters::TauTrackFlag flag=TauJetParameters::TauTrackFlag::classifiedCharged) const
@ RNNEleScore
RNN score for Ele rejection (not transformed)
Definition TauDefs.h:94
TauJet_v3 TauJet
Definition of the current "tau version".