PseudoJet.hh
Go to the documentation of this file.
1 
4 #pragma once
5 
6 #include <boost/operators.hpp>
7 
8 #include "fastjet/PseudoJet.hh"
9 #include "boca/units/Prefixes.hh"
10 
11 class TLorentzVector;
12 
13 namespace boca
14 {
15 
16 template<typename Value_>
18 
19 template<typename Value_>
20 class Vector3;
21 
22 template<typename Value_>
23 class Vector2;
24 
30 class PseudoJet : fastjet::PseudoJet
31  , boost::totally_ordered<PseudoJet>
32 {
33 
34 public:
35 
44  PseudoJet();
45 
49  PseudoJet(fastjet::PseudoJet const &pseudo_jet);
50 
54  PseudoJet(TLorentzVector const &vector);
55 
60 
64  PseudoJet(boca::Vector3<boca::Momentum> const& spatial, const boca::Energy &e);
65 
69  PseudoJet(Vector2<boca::Momentum> const& transverse, Momentum const& z, boca::Energy const& e);
70 
74  PseudoJet(Momentum const& x, Momentum const& y, Momentum const& z, boca::Energy const& e);
75 
76  PseudoJet(PseudoJet const&) = default;
77 
78  PseudoJet(PseudoJet &&) = default;
79 
81 
82  void Reset(PseudoJet const &pseudo_jet);
83 
89  void ScaleMomentum(double factor);
90 
99  fastjet::PseudoJet &FastJet();
100 
101  fastjet::PseudoJet const &FastJet() const;
102 
104 
109  virtual UserInfoBase const &Info() const;
110 
111  virtual UserInfoBase &Info();
112 
113  bool HasInfo();
114 
116 
125  Momentum Pt() const;
126 
130  boca::Energy Energy() const;
131 
135  boca::Energy E() const;
136 
140  boca::Energy Et2() const;
141 
145  boca::Energy Et() const;
146 
150  boca::Mass Mass() const;
151 
156 
160  Momentum Px() const;
161 
165  Momentum Py() const;
166 
170  Momentum Pz() const;
171 
175  MomentumSquare ModP2() const;
176 
180  Momentum ModP() const;
181 
183 
192  Angle Rap() const;
193 
197  Angle Phi() const;
198 
202  Angle DeltaPhiTo(PseudoJet const &jet) const;
203 
207  Angle DeltaRapTo(PseudoJet const &jet) const;
208 
212  Angle DeltaRTo(PseudoJet const &jet) const;
213 
215 
224  Vector2<Angle> DeltaTo(PseudoJet const &jet) const;
225 
230  Vector2<Angle> Angles(bool wrap_phi = false) const;
231 
235  Vector2<Angle> AnglesMinTo(PseudoJet const &jet) const;
236 
241 
246 
251 
253 
262  bool operator<(PseudoJet const &jet) const;
263 
267  bool operator==(const PseudoJet &pseudo_jet) const;
268 
272  PseudoJet &operator=(PseudoJet const&) & = default;
273 
277  PseudoJet &operator=(PseudoJet &&) & = default;
278 
282  friend std::ostream & operator<<(std::ostream & stream, PseudoJet const& jet);
283 
285 
286  virtual ~PseudoJet() {}
287 
288 protected:
289 
290 
291 };
292 
293 }
typename boost::units::multiply_typeof_helper< Mass, Mass >::type MassSquare
Definition: ElectronVolt.hh:75
Lorentz vector.
Definition: PseudoJet.hh:17
boca::Energy Energy() const
Energy.
Definition: PseudoJet.cpp:103
Vector2< Angle > AnglesMinTo(PseudoJet const &jet) const
Vector of rapidity and azimuth with minimal distance to jet.
Definition: PseudoJet.cpp:189
boca::Energy Et() const
transverse energy
Definition: PseudoJet.cpp:118
void Reset(PseudoJet const &pseudo_jet)
Definition: PseudoJet.cpp:224
boost::units::quantity< boost::units::si::plane_angle > Angle
Angle measured in radian.
Definition: Si.hh:35
Vector2< Angle > Angles(bool wrap_phi=false) const
Vector of rapidity and azimuth .
Definition: PseudoJet.cpp:184
Momentum Pz() const
Z component.
Definition: PseudoJet.cpp:133
boca::Vector3< Momentum > Spatial() const
Momentum three vector.
Definition: PseudoJet.cpp:78
boca::Energy Et2() const
Transverse energy square.
Definition: PseudoJet.cpp:113
friend std::ostream & operator<<(std::ostream &stream, PseudoJet const &jet)
Output stream operator.
Definition: PseudoJet.cpp:212
Angle DeltaRTo(PseudoJet const &jet) const
Difference to a jet.
Definition: PseudoJet.cpp:173
Wrapper for fastjet::PseudoJet adding BoCA related functions.
Definition: PseudoJet.hh:30
typename boost::units::multiply_typeof_helper< Momentum, Momentum >::type MomentumSquare
Definition: ElectronVolt.hh:76
PseudoJet()
Default Constructor.
Definition: PseudoJet.cpp:14
boca::Energy E() const
Energy.
Definition: PseudoJet.cpp:108
fastjet::PseudoJet & FastJet()
Fastjet PseudoJet.
Definition: PseudoJet.cpp:83
boca::MassSquare MassSquare() const
Mass square.
Definition: PseudoJet.cpp:123
virtual UserInfoBase const & Info() const
Definition: PseudoJet.cpp:61
bool operator==(const PseudoJet &pseudo_jet) const
Equality comparable.
Definition: PseudoJet.cpp:207
bool operator<(PseudoJet const &jet) const
Less than comparable.
Definition: PseudoJet.cpp:202
Vector2< Angle > DeltaTo(PseudoJet const &jet) const
Angular distance to a jet.
Definition: PseudoJet.cpp:179
Momentum Lorentz vector.
Definition: LorentzVector.hh:28
Angle Rap() const
Rapidity .
Definition: PseudoJet.cpp:158
boca::Mass Mass() const
Mass.
Definition: PseudoJet.cpp:98
Angle DeltaRapTo(PseudoJet const &jet) const
Difference to a jet.
Definition: PseudoJet.cpp:168
Angle Phi() const
Azimuth constrained to .
Definition: PseudoJet.cpp:153
Boosted Collider Analysis.
Definition: Analysis.hh:15
Momentum Px() const
X component.
Definition: PseudoJet.cpp:138
boca::LorentzVector< Momentum > LorentzVector() const
Momentum lorentz vector.
Definition: PseudoJet.cpp:73
Momentum ModP() const
the 3-vector modulus
Definition: PseudoJet.cpp:148
virtual ~PseudoJet()
Definition: PseudoJet.hh:286
Three dimensionial vector.
Definition: PseudoJet.hh:20
bool HasInfo()
Definition: PseudoJet.cpp:218
boost::units::quantity< electronvolt::Energy > Energy
Energy measured in electronvolt.
Definition: ElectronVolt.hh:56
void ScaleMomentum(double factor)
Rescale the jet momentum.
Definition: PseudoJet.cpp:56
Momentum Py() const
Y component.
Definition: PseudoJet.cpp:128
Vector2< Momentum > Transverse() const
Transverse momentum vector .
Definition: PseudoJet.cpp:197
Momentum Pt() const
Transverse Momentum.
Definition: PseudoJet.cpp:93
Angle DeltaPhiTo(PseudoJet const &jet) const
Difference to a jet constrained to .
Definition: PseudoJet.cpp:163
Energy Momentum
Momentum measured in electronvolt.
Definition: ElectronVolt.hh:68
MomentumSquare ModP2() const
squared 3-vector modulus
Definition: PseudoJet.cpp:143
PseudoJet & operator=(PseudoJet const &)&=default
Default copy constructor.
Two dimensional Vector.
Definition: PseudoJet.hh:23
Energy Mass
Mass measured in electronvolt.
Definition: ElectronVolt.hh:62