ATLAS Offline Software
Loading...
Searching...
No Matches
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
11namespace 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
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}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
Helper class to provide constant type-safe access to aux data.
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.
void fill(const xAOD::TauJet &tau, const xAOD::TruthParticle &, float weight)
ResolutionPlots(PlotBase *pParent, const std::string &sDir, std::string sTauJetContainerName)
virtual double phi() const
The azimuthal angle ( ) of the particle.
float charge() const
virtual double pt() const
The transverse momentum ( ) of the particle.
virtual double m() const
The invariant mass of the particle.
virtual double eta() const
The pseudorapidity ( ) of the particle.
STL namespace.
TauJet_v3 TauJet
Definition of the current "tau version".
TruthParticle_v1 TruthParticle
Typedef to implementation.