ATLAS Offline Software
PhysicsAnalysis/TauID/TauDQA/src/ResolutionPlots.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 #include <utility>
6 #include "ResolutionPlots.h"
8 #include "TLorentzVector.h"
10 
11 namespace Tau{
12 
13  ResolutionPlots::ResolutionPlots(PlotBase *pParent, const std::string& sDir, std::string sTauJetContainerName):
14  PlotBase(pParent, sDir),
15  m_sTauJetContainerName(std::move(sTauJetContainerName))
16  {
17  }
18 
20  {
21  }
22 
23 
25  m_ptResolution = Book1D("ptResolution",m_sTauJetContainerName + " pt Resolution; pt(TauReco)/pt(visTauTruth); # Part",100,0.5,1.5);
26  m_etaResolution = Book1D("etaResolution",m_sTauJetContainerName + " eta Resolution; eta(TauReco) - eta(visTauTruth); # Part",100,-0.1,0.1);
27  m_phiResolution = Book1D("phiResolution",m_sTauJetContainerName + " phi Resolution; phi(TauReco) - phi(visTauTruth); # Part",100,-0.1,0.1);
28  m_chargeResolution = Book1D("chargeResolution",m_sTauJetContainerName + " charge Resolution; charge(TauReco) - charge(TauTruth); # Part",12,-6,6);
29  }
30 
31  void ResolutionPlots::fill(const xAOD::TauJet& tau, const xAOD::TruthParticle& truthtau, float weight) {
32 
33  static const SG::ConstAccessor<double> acc_ptvis("pt_vis");
34  static const SG::ConstAccessor<double> acc_etavis("eta_vis");
35  static const SG::ConstAccessor<double> acc_phivis("phi_vis");
36  static const SG::ConstAccessor<double> acc_mvis("m_vis");
37  static const SG::ConstAccessor<int> acc_pdgID("pdgId");
38  float ptratio = -999;
39  float pt = acc_ptvis(truthtau);
40  if(pt>0.) ptratio = tau.pt()/pt;
41  float eta = acc_etavis(truthtau);
42  float phi = acc_phivis(truthtau);
43  float m = acc_mvis(truthtau);
44  m_ptResolution->Fill(ptratio, weight);
45  m_etaResolution->Fill(tau.eta() - eta, weight);
46  TLorentzVector LV_TauJet(0,0,0,0);
47  TLorentzVector LV_TruthTau(0,0,0,0);
48  LV_TauJet.SetPtEtaPhiM(tau.pt(), tau.eta(), tau.phi(), tau.m());
49  LV_TruthTau.SetPtEtaPhiM(pt, eta, phi, m);
50  m_phiResolution->Fill(LV_TauJet.DeltaPhi(LV_TruthTau), weight);
51 
52  int pdgID = acc_pdgID(truthtau);
53  float truth_charge = 0.;
54  if(pdgID > 0) truth_charge = -1.;
55  else if(pdgID < 0) truth_charge = 1.;
56  float charge_diff = tau.charge() - truth_charge;
57 
58  m_chargeResolution->Fill(charge_diff, weight);
59 
60  }
61 
62 }
phi
Scalar phi() const
phi method
Definition: AmgMatrixBasePlugin.h:67
PlotBase
Definition: PlotBase.h:34
xAOD::TauJet_v3::eta
virtual double eta() const
The pseudorapidity ( ) of the particle.
Tau::ResolutionPlots::m_sTauJetContainerName
std::string m_sTauJetContainerName
Definition: PhysicsAnalysis/TauID/TauDQA/src/ResolutionPlots.h:29
ResolutionPlots.h
eta
Scalar eta() const
pseudorapidity method
Definition: AmgMatrixBasePlugin.h:83
python.copyTCTOutput.sDir
sDir
Definition: copyTCTOutput.py:57
Tau::ResolutionPlots::ResolutionPlots
ResolutionPlots(PlotBase *pParent, const std::string &sDir, std::string sTauJetContainerName)
Definition: PhysicsAnalysis/TauID/TauDQA/src/ResolutionPlots.cxx:13
xAOD::TauJet_v3::m
virtual double m() const
The invariant mass of the particle.
Tau::ResolutionPlots::~ResolutionPlots
virtual ~ResolutionPlots()
Definition: PhysicsAnalysis/TauID/TauDQA/src/ResolutionPlots.cxx:19
Tau::ResolutionPlots::fill
void fill(const xAOD::TauJet &tau, const xAOD::TruthParticle &, float weight)
Definition: PhysicsAnalysis/TauID/TauDQA/src/ResolutionPlots.cxx:31
test_pyathena.pt
pt
Definition: test_pyathena.py:11
Tau
Definition: TauDQA/src/CorePlots.cxx:8
Tau::ResolutionPlots::initializePlots
void initializePlots()
Definition: PhysicsAnalysis/TauID/TauDQA/src/ResolutionPlots.cxx:24
SG::ConstAccessor< double >
dqt_zlumi_pandas.weight
int weight
Definition: dqt_zlumi_pandas.py:190
PlotBase::Book1D
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
xAOD::TauJet_v3::pt
virtual double pt() const
The transverse momentum ( ) of the particle.
xAOD::TauJet_v3::phi
virtual double phi() const
The azimuthal angle ( ) of the particle.
ITauTruthMatchingTool.h
Tau::ResolutionPlots::m_phiResolution
TH1 * m_phiResolution
Definition: PhysicsAnalysis/TauID/TauDQA/src/ResolutionPlots.h:24
xAOD::TruthParticle_v1
Class describing a truth particle in the MC record.
Definition: TruthParticle_v1.h:37
xAOD::TauJet_v3
Class describing a tau jet.
Definition: TauJet_v3.h:41
Tau::ResolutionPlots::m_etaResolution
TH1 * m_etaResolution
Definition: PhysicsAnalysis/TauID/TauDQA/src/ResolutionPlots.h:23
xAOD::TauJet_v3::charge
float charge() const
Tau::ResolutionPlots::m_ptResolution
TH1 * m_ptResolution
Definition: PhysicsAnalysis/TauID/TauDQA/src/ResolutionPlots.h:22
Tau::ResolutionPlots::m_chargeResolution
TH1 * m_chargeResolution
Definition: PhysicsAnalysis/TauID/TauDQA/src/ResolutionPlots.h:25
ConstAccessor.h
Helper class to provide constant type-safe access to aux data.
python.SystemOfUnits.m
float m
Definition: SystemOfUnits.py:106