ATLAS Offline Software
Loading...
Searching...
No Matches
egammaConversionFillerTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration
3*/
4
5// $Id$
12
13
15#include "xAODTracking/Vertex.h"
17#include "CLHEP/Vector/LorentzVector.h"
18#include <cmath>
19
20
21using CLHEP::HepLorentzVector;
22
23
24namespace D3PD {
25
26
34 (const std::string& type,
35 const std::string& name,
36 const IInterface* parent)
37 : BlockFillerTool<xAOD::Photon> (type, name, parent)
38{
39 book().ignore(); // Avoid coverity warnings.
40}
41
42
47{
48 CHECK( addVariable ("hasconv", m_hasconv) );
49 CHECK( addVariable ("convvtxx", m_convvtxx) );
50 CHECK( addVariable ("convvtxy", m_convvtxy) );
51 CHECK( addVariable ("convvtxz", m_convvtxz) );
52 CHECK( addVariable ("Rconv", m_Rconv) );
53 CHECK( addVariable ("zconv", m_zconv) );
54 CHECK( addVariable ("convvtxchi2", m_convvtxchi2) );
55 CHECK( addVariable ("pt1conv", m_pt1conv) );
56 CHECK( addVariable ("convtrk1nBLHits", m_convtrk1nBLHits) );
57 CHECK( addVariable ("convtrk1nPixHits", m_convtrk1nPixHits) );
58 CHECK( addVariable ("convtrk1nSCTHits", m_convtrk1nSCTHits) );
59 CHECK( addVariable ("convtrk1nTRTHits", m_convtrk1nTRTHits) );
60 CHECK( addVariable ("pt2conv", m_pt2conv) );
61 CHECK( addVariable ("convtrk2nBLHits", m_convtrk2nBLHits) );
62 CHECK( addVariable ("convtrk2nPixHits", m_convtrk2nPixHits) );
63 CHECK( addVariable ("convtrk2nSCTHits", m_convtrk2nSCTHits) );
64 CHECK( addVariable ("convtrk2nTRTHits", m_convtrk2nTRTHits) );
65 CHECK( addVariable ("ptconv", m_ptconv) );
66 CHECK( addVariable ("pzconv", m_pzconv) );
67
68 return StatusCode::SUCCESS;
69}
70
71
72namespace {
73
74
76const xAOD::TrackParticle* gettp (const xAOD::Vertex* conv, unsigned int n)
77{
78 if (n >= conv->nTrackParticles()) return 0;
79 return conv->trackParticle(n);
80}
81
82
83}
84
85
97{
98 const xAOD::Vertex* conv = p.vertex();
99 if (conv) {
100 *m_hasconv = true;
101
102 *m_convvtxx = conv->x();
103 *m_convvtxy = conv->y();
104 *m_convvtxz = conv->z();
105 *m_Rconv = static_cast<float> (hypot (conv->x(), conv->y()));
106 *m_zconv = static_cast<float> (conv->z());
107 *m_convvtxchi2 = static_cast<float> (conv->chiSquared());
108
109 const xAOD::TrackParticle* tp1 = gettp (conv, 0);
110 const xAOD::TrackParticle* tp2 = gettp (conv, 1);
111 TLorentzVector psum;
112 if (tp1) {
113 psum += tp1->p4();
114 *m_pt1conv = tp1->pt();
119 }
120
121
122 if (tp2) {
123 psum += tp2->p4();
124 *m_pt2conv = tp2->pt();
129 }
130
131 *m_ptconv = psum.Pt();
132 *m_pzconv = psum.Pz();
133 }
134
135 return StatusCode::SUCCESS;
136}
137
138
139} // namespace D3PD
Helpers for checking error return status codes and reporting errors.
#define CHECK(...)
Evaluate an expression and check for errors.
virtual StatusCode addVariable(const std::string &name, const std::type_info &ti, void *&ptr, const std::string &docstring="", const void *defval=0)
Type-safe wrapper for block filler tools.
uint8_t * m_convtrk1nTRTHits
Variable: conversion track 1 number of TRT hits.
bool * m_hasconv
Variable: is there a conversion?
uint8_t * m_convtrk1nSCTHits
Variable: conversion track 1 number of SCT hits.
uint8_t * m_convtrk1nBLHits
Variable: conversion track 1 number of B layer hits.
float * m_ptconv
Variable: pt of conversion tracks 1+2.
uint8_t * m_convtrk2nTRTHits
Variable: conversion track 2 number of TRT hits.
egammaConversionFillerTool(const std::string &type, const std::string &name, const IInterface *parent)
Standard Gaudi tool constructor.
float * m_convvtxy
Variable: conversion vertex y.
float * m_Rconv
Variable: radius of conversion.
virtual StatusCode book() final
Book variables for this block.
float * m_convvtxx
Variable: conversion vertex x.
float * m_zconv
Variable: z of conversion.
uint8_t * m_convtrk2nPixHits
Variable: conversion track 2 number of pixel hits.
uint8_t * m_convtrk1nPixHits
Variable: conversion track 1 number of pixel hits.
float * m_pt2conv
Variable: pt of track2 of conversion.
virtual StatusCode fill(const xAOD::Photon &p) override
Fill one block — type-safe version.
float * m_pt1conv
Variable: pt of track1 of conversion.
float * m_convvtxz
Variable: conversion vertex z.
float * m_convvtxchi2
Variable: conversion vertex chi2.
uint8_t * m_convtrk2nBLHits
Variable: conversion track 2 number of B layer hits.
uint8_t * m_convtrk2nSCTHits
Variable: conversion track 2 number of SCT hits.
float * m_pzconv
Variable: pz of conversion tracks 1+2.
Class describing an photon.
virtual FourMom_t p4() const override final
The full 4-momentum of the particle.
bool summaryValue(uint8_t &value, const SummaryType &information) const
Accessor for TrackSummary values.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
Block filler tool for photon conversion information.
Block filler tool for noisy FEB information.
ICaloAffectedTool is abstract interface for tools checking if 4 mom is in calo affected region.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Vertex_v1 Vertex
Define the latest version of the vertex class.
Photon_v1 Photon
Definition of the current "egamma version".
@ numberOfTRTHits
number of TRT hits [unit8_t].
@ numberOfSCTHits
number of hits in SCT [unit8_t].
@ numberOfInnermostPixelLayerHits
these are the hits in the 0th pixel barrel layer
@ numberOfPixelHits
these are the pixel hits, including the b-layer [unit8_t].