ATLAS Offline Software
Loading...
Searching...
No Matches
FlowElement_v1.h
Go to the documentation of this file.
1// this files is -*- c++ -*-
2/*
3 Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
4*/
5#ifndef XAODPFLOW_VERSIONS_FLOWELEMENT_V1_H
6#define XAODPFLOW_VERSIONS_FLOWELEMENT_V1_H
7
8// Core include(s):
9#include "AthLinks/ElementLink.h"
10
11// xAOD include(s):
12#include "xAODBase/IParticle.h"
14
15namespace xAOD {
16
25 class FlowElement_v1 : public IParticle {
26 public:
27
29
30 typedef unsigned long signal_t; // 32-bit minimum
31 typedef short vertex_t;
33
35 enum SignalType {
36 // global characteristics
37 Neutral = 0x1000,
38 Charged = 0x2000,
39 Combined = 0x4000, //needed??
40
41 // detector level signals and flow objects
42 CaloCluster = Neutral | 0x0100,
43 Track = Charged | 0x0200,
44 Muon = Charged | 0x0400,
45 PFlow = 0x0010,
48
49 // higher level flow objects
50 TCC = 0x0020,
53 UFO = 0x0001,
56
57 // unknown
58 Unknown = 0x0000
59 };
60
63 Undefined = 0x00, HardScatter = 0x10, Pileup = 0x20, PileupSideBand = 0x21
64 };
65
66 // *************************************************
69 virtual double pt() const override;
70 virtual double eta() const override;
71 virtual double phi() const override;
72 virtual double m() const override;
73 virtual double e() const override;
74 virtual double rapidity() const override;
75 virtual FourMom_t p4() const override;
76 virtual Type::ObjectType type() const override ;
77
78 void setP4(float pt, float eta, float phi, float m) ;
79 void setP4(const FourMom_t& p4);
81
82 // *************************************************
88
89 // *************************************************
92 bool isMatchedToPV(MatchedPVType vxtype=HardScatter) const;
93
97
98 // *************************************************
101 bool isCharged() const;
102 float charge() const;
103
104 void setCharge(float c);
106
107 // *************************************************
108
112 std::vector<const xAOD::IParticle*> chargedObjects() const ;
113 std::vector<std::pair<const xAOD::IParticle*,float> > chargedObjectsAndWeights() const ;
114
115 std::size_t nChargedObjects() const ;
116 const xAOD::IParticle* chargedObject( std::size_t i ) const;
117 std::pair< const xAOD::IParticle*, float > chargedObjectAndWeight( std::size_t i ) const;
119
120
122
123 const std::vector<ElementLink<IParticleContainer>> & chargedObjectLinks() const ;
124 const std::vector<float> & chargedObjectWeights() const ;
126 void setChargedObjectLinks(const std::vector<ElementLink<IParticleContainer>> & elV, const std::vector<float> & wV);
128
129 // *************************************************
133 std::vector<const xAOD::IParticle*> otherObjects() const ;
134 std::vector<std::pair<const xAOD::IParticle*,float> > otherObjectsAndWeights() const ;
135
136 std::size_t nOtherObjects() const ;
137 const xAOD::IParticle* otherObject( std::size_t i ) const;
138 std::pair< const xAOD::IParticle*, float > otherObjectAndWeight( std::size_t i ) const;
140
143 const std::vector<ElementLink<IParticleContainer>>& otherObjectLinks() const ;
144 const std::vector<float>& otherObjectWeights() const ;
145
147 void setOtherObjectLinks(const std::vector<ElementLink<IParticleContainer>> & elV, const std::vector<float> & wV);
149
150 };
151
152}
153// Declare IParticle as a base class of FlowElement_v1:
156
157
158#endif
An STL vector of pointers that by default owns its pointed-to elements.
#define DATAVECTOR_BASE(T, BASE)
Declare base class info to DataVector.
Definition DataVector.h:649
A detector object made of other lower level object(s)
std::size_t nChargedObjects() const
std::size_t nOtherObjects() const
std::vector< const xAOD::IParticle * > chargedObjects() const
virtual double rapidity() const override
The true rapidity (y) of the particle.
std::vector< const xAOD::IParticle * > otherObjects() const
const std::vector< float > & chargedObjectWeights() const
void setVertexType(vertex_t t)
vertex_t vertexType() const
virtual double pt() const override
void setP4(float pt, float eta, float phi, float m)
virtual double m() const override
The invariant mass of the particle.
bool isMatchedToPV(MatchedPVType vxtype=HardScatter) const
const std::vector< ElementLink< IParticleContainer > > & otherObjectLinks() const
virtual double phi() const override
The azimuthal angle ( ) of the particle.
void setChargedObjectLinks(const std::vector< ElementLink< IParticleContainer > > &elV, const std::vector< float > &wV)
MatchedPVType
Enum to encode high-level information on the vertex associated to this FlowElement.
void setCharge(float c)
const std::vector< ElementLink< IParticleContainer > > & chargedObjectLinks() const
Access to the EL.
virtual double eta() const override
The pseudorapidity ( ) of the particle.
void setOtherObjectLinks(const std::vector< ElementLink< IParticleContainer > > &elV, const std::vector< float > &wV)
signal_t signalType() const
std::vector< std::pair< const xAOD::IParticle *, float > > otherObjectsAndWeights() const
void setChargedObjectLinks(const std::vector< ElementLink< IParticleContainer > > &elV)
std::pair< const xAOD::IParticle *, float > chargedObjectAndWeight(std::size_t i) const
void setSignalType(signal_t t)
TLorentzVector FourMom_t
Definition of the 4-momentum type.
const std::vector< float > & otherObjectWeights() const
SignalType
Enum to encode the nature of the object this FlowElement represents.
const xAOD::IParticle * chargedObject(std::size_t i) const
float charge() const
unsigned long signal_t
virtual Type::ObjectType type() const override
The type of the object as a simple enumeration.
void setOtherObjectLinks(const std::vector< ElementLink< IParticleContainer > > &elV)
virtual double e() const override
The total energy of the particle.
const xAOD::IParticle * otherObject(std::size_t i) const
std::vector< std::pair< const xAOD::IParticle *, float > > chargedObjectsAndWeights() const
virtual FourMom_t p4() const override
The full 4-momentum of the particle.
std::pair< const xAOD::IParticle *, float > otherObjectAndWeight(std::size_t i) const
Class providing the definition of the 4-vector interface.
IParticle()=default
TLorentzVector FourMom_t
Definition of the 4-momentum type.
NRpcCablingAlg reads raw condition data and writes derived condition data to the condition store.
Definition PFlow.py:1
ICaloAffectedTool is abstract interface for tools checking if 4 mom is in calo affected region.