ATLAS Offline Software
Loading...
Searching...
No Matches
MomentumPullPlots.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration
3*/
4
9
10using namespace xAOD::P4Helpers;
11namespace Muon{
12
14
15 dpt_idme = Book1D("dpt_idme", "idme_dpt;pt_{id} - pt_{me} [GeV];Entries / 0.2 GeV", 100,-100.,100);
16 ddpt_idme = Book1D("ddpt_idme", "idme_ddpt;(pt_{id} - pt_{me})/pt_{id};Entries", 300,-0.5,0.5);
17 dphi_idme = Book1D("dphi_idme", "idme_dphi;#phi_{id} - #phi_{me};Entries", 100,-0.1,0.1);
18 deta_idme = Book1D("deta_idme", "idme_deta;#eta_{id} - #eta_{me};Entries", 100,-0.1,0.1);
19
20 dpt_cbme = Book1D("dpt_cbme", "cbme_dpt;pt_{cb} - pt_{me} [GeV];Entries / 0.2 GeV", 100,-100.,100);
21 ddpt_cbme = Book1D("ddpt_cbme", "cbme_ddpt;(pt_{cb} - pt_{me})/pt_{cb};Entries", 300,-0.5,0.5);
22 dphi_cbme = Book1D("dphi_cbme", "cbme_dphi;#phi_{cb} - #phi_{me};Entries", 100,-0.1,0.1);
23 deta_cbme = Book1D("deta_cbme", "cbme_deta;#eta_{cb} - #eta_{me};Entries", 100,-0.1,0.1);
24
25 pt_cbme = Book2D("pt_cbme", "pt_cbme; pt_{cb} [GeV]; pt_{me} [GeV]", 100, 0.0, 100.0, 100, 0.0, 100.0);
26 pt_cbid = Book2D("pt_cbid", "pt_cbid; pt_{cb} [GeV]; pt_{id} [GeV]", 100, 0.0, 100.0, 100, 0.0, 100.0);
27 pt_meid = Book2D("pt_meid", "pt_meid; pt_{me} [GeV]; pt_{id} [GeV]", 100, 0.0, 100.0, 100, 0.0, 100.0);
28}
29
30 void MomentumPullPlots::fill(const xAOD::Muon& mu, float weight) {
31 const xAOD::TrackParticle* cb = mu.trackParticle(xAOD::Muon::CombinedTrackParticle);
32 const xAOD::TrackParticle* id = mu.trackParticle(xAOD::Muon::InnerDetectorTrackParticle);
33
35
36 //check correct numbering in Muon.h -> OK!
37 const xAOD::TrackParticle* me = mu.trackParticle(xAOD::Muon::TrackParticleType::ExtrapolatedMuonSpectrometerTrackParticle);
38 if (!me) me = mu.trackParticle(xAOD::Muon::MuonSpectrometerTrackParticle);
39
40 if (cb && me){
41 dpt_cbme->Fill( (cb->pt() - me->pt())*0.001, weight);
42 ddpt_cbme->Fill( (cb->pt() - me->pt())/cb->pt(), weight);
43 dphi_cbme->Fill( deltaPhi(cb->phi(), me->phi()),weight);
44 deta_cbme->Fill( cb->eta() - me->eta(), weight);
45 pt_cbme->Fill( cb->pt()*0.001, me->pt()*0.001 , weight);
46 }
47
48 if (id && me){
49 dpt_idme->Fill( (id->pt() - me->pt())*0.001, weight);
50 ddpt_idme->Fill( (id->pt() - me->pt())/id->pt(), weight);
51 dphi_idme->Fill( deltaPhi(id->phi() , me->phi()), weight);
52 deta_idme->Fill( id->eta() - me->eta(), weight);
53 pt_meid->Fill( me->pt()*0.001, id->pt()*0.001 , weight);
54 }
55
56 if(id && cb)
57 pt_cbid->Fill( cb->pt()*0.001, id->pt()*0.001, weight );
58}
59
60}
Scalar deltaPhi(const MatrixBase< Derived > &vec) const
void fill(const xAOD::Muon &mu, float weight=1.0)
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
TH2F * Book2D(const std::string &name, const std::string &labels, int nBinsX, float startX, float endX, int nBinsY, float startY, float endY, bool prependDir=true)
Book a TH2F histogram.
Definition PlotBase.cxx:123
virtual double phi() const override final
The azimuthal angle ( ) of the particle (has range to .)
virtual double pt() const override final
The transverse momentum ( ) of the particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
NRpcCablingAlg reads raw condition data and writes derived condition data to the condition store.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Muon_v1 Muon
Reference the current persistent version: