ATLAS Offline Software
Loading...
Searching...
No Matches
JetVariable.h
Go to the documentation of this file.
1// -*- c++ -*-
2
3/*
4 Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
5*/
6
7#ifndef JETMONITORING_JETVARIABLE_H
8#define JETMONITORING_JETVARIABLE_H
9#include <vector>
10
11#include "xAODJet/Jet.h"
12
26
27
28
29namespace JetVar {
30 template<typename T>
32
33
37 struct VectorValue {
39 virtual ~VectorWrapper(){}
40 virtual float at(int) const =0;
41 virtual size_t size() const =0;
42 };
43 VectorValue(VectorWrapper * v=nullptr , float s=1.) : m_v(v), m_scale(s){}
44 float operator[](int i) const {return m_v->at(i)*m_scale;}
45 size_t size()const {return m_v->size();}
46 std::unique_ptr<VectorWrapper> m_v;
47 float m_scale=1;
48 };
49
50
51 class Variable {
52 public:
53
54 Variable(const std::string & name) : m_name(name) {}
55 virtual ~Variable() = default;
56 virtual float value(const xAOD::Jet &) const = 0;
57
59 virtual bool isVector() const {return false;}
60
62 virtual VectorValue vector(const xAOD::Jet &) const {return VectorValue();}
63
64
65 virtual std::string name() const {return m_name;}
66
67 float scale() const {return m_scale;}
68 void setScale(float s) { m_scale=s;}
69
70
71
74 static std::unique_ptr<Variable> create(const std::string & name, const std::string &type="float", int index=-1);
75
76 std::string m_name;
77 float m_scale = 1;
78 int m_index = -1;
79 };
80
81
82
83
84 // ********************************************************
85 // Concrete implementations of Variable
86
87
89 template<typename T>
90 struct VariableAtt : public Variable {
91 VariableAtt(const std::string & name) : Variable(name), m_acc(name) {}
92 virtual float value(const xAOD::Jet & j) const {
93 if ( m_acc.isAvailable( j ) ) return m_acc(j)*m_scale;
94 else return -999.;
95 }
97 };
98
99
100
102 template<typename T>
103 struct VariableAtt<std::vector<T> > : public Variable {
104 typedef typename std::vector<T> vect_t;
106 VectorWrapperT(const vect_t *v): m_v(v){};
107 float at(int i) const {return (*m_v)[i];}
108 size_t size() const {return m_v->size();}
109 const vect_t * m_v;
110 };
111
112 VariableAtt(const std::string & name, int index) : Variable(name), m_acc(name) {
114 }
115
116 virtual std::string name() const {
117 if(isVector() )return m_name;
118 return m_name+std::to_string(m_index);
119 }
120
121 // returns false if the index is valid : in this case it is a simple variable
122 virtual bool isVector() const {return m_index==-1;}
123
124 // use only if the index is valid
125 virtual float value(const xAOD::Jet & j) const {
126 if ( m_acc.isAvailable( j ) && m_index >= 0 && m_acc(j).size() > (unsigned) m_index ) {
127 return m_acc(j)[m_index]*m_scale;
128 }
129 else return -999.;
130 }
131
132 virtual VectorValue vector(const xAOD::Jet &j) const {
133 if ( m_acc.isAvailable( j ) ) {
134 VectorValue v( new VectorWrapperT(&m_acc(j)) , m_scale ) ;
135 return v;
136 }
137 else {
139 return junk;
140 }
141 }
142
143
145 // this is the value used returned when a vector attribute is not accessible
146 const static inline vect_t s_dummy{-999};
147 };
148
149
150
151
152 // *******************************************
153 // The classes below represent variables not stored as attribute in the Jet EDM
154 //
155
156 struct EVar : public Variable {
157 using Variable::Variable;
158 virtual float value(const xAOD::Jet & j) const { return j.e()*m_scale;}
159 };
160
161 struct PzVar : public Variable {
162 using Variable::Variable;
163 virtual float value(const xAOD::Jet & j) const { return j.pz()*m_scale;}
164 };
165
166 struct NconstitVar : public Variable {
167 using Variable::Variable;
168 virtual float value(const xAOD::Jet & j) const { return j.numConstituents();}
169 };
170
171 struct Rapidity : public Variable {
172 using Variable::Variable;
173 virtual float value(const xAOD::Jet & j) const { return j.rapidity();}
174 };
175
176 struct AbsEtaVar : public Variable {
177 using Variable::Variable;
178 virtual float value(const xAOD::Jet & j) const { return fabs(j.eta());}
179 };
180
181 struct EtVar : public Variable {
182 using Variable::Variable;
183 virtual float value(const xAOD::Jet & j) const { return j.p4().Et()*m_scale;}
184 };
185
186 struct FChargedVar : public Variable {
187 using Variable::Variable;
188 virtual float value(const xAOD::Jet & j) const {
189 bool status = false;
190 float constScalePt = 0.;
191 std::vector<float> SumPtChargedPFOPt500;
192 status = j.getAttribute<float>("JetConstitScaleMomentum_pt", constScalePt ); // Jet pT at the constituent scale
193 if (!status) return 0;
194 status = j.getAttribute<std::vector<float> >("SumPtChargedPFOPt500", SumPtChargedPFOPt500 ); //Vector over all vertices in the event, each element contains the sum pT of all charged PFO with a pT > 0.5 GeV associated to the vertex.
195 if (!status) return 0;
196 return SumPtChargedPFOPt500.at(0)/=constScalePt; //definition of "fCharge", index 0 points to the primary vertex
197 }
198 };
199
200 struct EM3FracVar : public Variable {
201 using Variable::Variable;
202 virtual float value(const xAOD::Jet & j) const {
203 float constitScaleEnergy = 0.;
204 std::vector<float> samplingFrac;
205 xAOD::JetFourMom_t fourVec;
206 bool status = false;
207
208 status = j.getAttribute<xAOD::JetFourMom_t>( "JetConstitScaleMomentum", fourVec ); // Jet four-momentum at constituent scale
209 if( status ) constitScaleEnergy = fourVec.E() * m_scale ;
210 else return 0.;
211 status = j.getAttribute<std::vector<float> >("EnergyPerSampling", samplingFrac ); //EnergyPerSampling is a vector of size 24; element i refers to the energy deposited in calo sampling i, see https://twiki.cern.ch/twiki/bin/viewauth/AtlasProtected/Run2JetMoments#Sampling_layers
212 if( status ) return (samplingFrac[3]+samplingFrac[7])/constitScaleEnergy; //3 is 'EMB3' in the LAr barrel, 7 is 'EME3' in the LAr EM endcap
213 else return 0.;
214 }
215 };
216
217 struct Tile0FracVar : public Variable {
218 using Variable::Variable;
219 virtual float value(const xAOD::Jet & j) const {
220 float constitScaleEnergy = 0.;
221 std::vector<float> samplingFrac;
222 xAOD::JetFourMom_t fourVec;
223 bool status = false;
224
225 status = j.getAttribute<xAOD::JetFourMom_t>( "JetConstitScaleMomentum", fourVec ); // Jet four-momentum at constituent scale
226 if( status ) constitScaleEnergy = fourVec.E() * m_scale ;
227 else return 0.;
228 status = j.getAttribute<std::vector<float> >("EnergyPerSampling", samplingFrac ); // refer to EM3FracVar above
229 if( status ) return (samplingFrac[12]+samplingFrac[18])/constitScaleEnergy; //12 is 'TileBar0' in the Tile barrel, 18 is 'TileExt0' in the Tile extended barrel
230 else return 0.;
231 }
232 };
233
234}
235
236#endif
float scale() const
Definition JetVariable.h:67
virtual std::string name() const
Definition JetVariable.h:65
virtual ~Variable()=default
virtual VectorValue vector(const xAOD::Jet &) const
return a helper object allowing to iterate over the underlying vector. !!! use only if isVector()==tr...
Definition JetVariable.h:62
Variable(const std::string &name)
Definition JetVariable.h:54
void setScale(float s)
Definition JetVariable.h:68
std::string m_name
Definition JetVariable.h:76
virtual float value(const xAOD::Jet &) const =0
static std::unique_ptr< Variable > create(const std::string &name, const std::string &type="float", int index=-1)
create and return a new Variable of a given name & type.
virtual bool isVector() const
return true if the underlying variable is of type vector<X>
Definition JetVariable.h:59
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
size_t numConstituents() const
Number of constituents in this jets (this is valid even when reading a file where the constituents ha...
Definition Jet_v1.cxx:153
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition Jet_v1.cxx:71
virtual double rapidity() const
The true rapidity (y) of the particle.
Definition Jet_v1.cxx:67
float pz() const
The z-component of the jet's momentum.
Definition Jet_v1.cxx:99
bool getAttribute(AttributeID type, T &value) const
Retrieve attribute moment by enum.
virtual double eta() const
The pseudorapidity ( ) of the particle.
Definition Jet_v1.cxx:49
virtual double e() const
The total energy of the particle.
Definition Jet_v1.cxx:63
SG::AuxElement::Accessor< T > Accessor
Definition JetVariable.h:31
Definition index.py:1
STL namespace.
Jet_v1 Jet
Definition of the current "jet version".
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiM4D< double > > JetFourMom_t
Base 4 Momentum type for Jet.
Definition JetTypes.h:17
Variable(const std::string &name)
Definition JetVariable.h:54
virtual float value(const xAOD::Jet &j) const
Variable(const std::string &name)
Definition JetVariable.h:54
virtual float value(const xAOD::Jet &j) const
Variable(const std::string &name)
Definition JetVariable.h:54
virtual float value(const xAOD::Jet &j) const
Variable(const std::string &name)
Definition JetVariable.h:54
virtual float value(const xAOD::Jet &j) const
virtual float value(const xAOD::Jet &j) const
Variable(const std::string &name)
Definition JetVariable.h:54
Variable(const std::string &name)
Definition JetVariable.h:54
virtual float value(const xAOD::Jet &j) const
Variable(const std::string &name)
Definition JetVariable.h:54
virtual float value(const xAOD::Jet &j) const
virtual float value(const xAOD::Jet &j) const
Variable(const std::string &name)
Definition JetVariable.h:54
Variable(const std::string &name)
Definition JetVariable.h:54
virtual float value(const xAOD::Jet &j) const
virtual VectorValue vector(const xAOD::Jet &j) const
return a helper object allowing to iterate over the underlying vector. !!! use only if isVector()==tr...
VariableAtt(const std::string &name, int index)
virtual std::string name() const
virtual bool isVector() const
return true if the underlying variable is of type vector<X>
virtual float value(const xAOD::Jet &j) const
virtual float value(const xAOD::Jet &j) const
Definition JetVariable.h:92
VariableAtt(const std::string &name)
Definition JetVariable.h:91
Accessor< T > m_acc
Definition JetVariable.h:96
virtual size_t size() const =0
virtual float at(int) const =0
VectorValue is a helper class to access any jet variable of type vector<X> It is implemented this way...
Definition JetVariable.h:37
VectorValue(VectorWrapper *v=nullptr, float s=1.)
Definition JetVariable.h:43
std::unique_ptr< VectorWrapper > m_v
Definition JetVariable.h:46
size_t size() const
Definition JetVariable.h:45
float operator[](int i) const
Definition JetVariable.h:44