Multiplet.hh
Go to the documentation of this file.
1 
4 #pragma once
5 
11 
12 namespace boca
13 {
14 
21 class Multiplet : public Identification
22 {
23 
24  template<typename Multiplet_>
25  using NotJet = typename std::enable_if < !std::is_same<Multiplet_, boca::Jet>::value && !std::is_same<Multiplet_, boca::PseudoJet>::value && !std::is_same<Multiplet_, boca::Particle>::value >::type;
26 
27 public:
28 
29  void SetClosestLepton(std::vector<boca::Lepton> const &leptons);
30 
31  ClosestLepton Lepton() const;
32 
33  boca::Jet Jet() const;
34 
35  boca::Mass Mass() const;
36 
37  Momentum Pt() const;
38 
47  Angle Rap() const;
48 
52  Angle Phi() const;
53 
57  template <typename Multiplet_, typename = NotJet<Multiplet_>>
58  Angle DeltaPhiTo(Multiplet_ const &multiplet) const
59  {
60  return Jet().DeltaPhiTo(multiplet.Jet());
61  }
62  Angle DeltaPhiTo(PseudoJet const &jet) const;
63 
67  template <typename Multiplet_, typename = NotJet<Multiplet_>>
68  Angle DeltaRapTo(Multiplet_ const &multiplet) const
69  {
70  return Jet().DeltaRapTo(multiplet.Jet());
71  }
72  Angle DeltaRapTo(PseudoJet const &jet) const;
73 
77  template <typename Multiplet_, typename = NotJet<Multiplet_>>
78  Angle DeltaRTo(Multiplet_ const &multiplet) const
79  {
80  return Jet().DeltaRTo(multiplet.Jet());
81  }
82  Angle DeltaRTo(PseudoJet const &jet) const;
83 
87  template <typename Multiplet_, typename = NotJet<Multiplet_>>
88  Vector2<Angle> DeltaTo(Multiplet_ const &multiplet) const
89  {
90  return Jet().DeltaTo(multiplet.Jet());
91  }
92  Vector2<Angle> DeltaTo(PseudoJet const &jet) const;
93 
98  Vector2<Angle> Angles(bool wrap_phi = false) const;
99 
103  template <typename Multiplet_/*, typename = NotJet<Multiplet_>*/>
104  Vector2<Angle> AnglesMinTo(Multiplet_ const &multiplet) const
105  {
106  return Jet().AnglesMinTo(multiplet.Jet());
107  }
108  Vector2<Angle> AnglesMinTo(PseudoJet const& jet) const;
110 
120 
124  bool HasConstituents() const;
125 
129  std::vector<boca::Jet> Constituents() const;
130 
136 
137  void SetExtraInfo(double extra_info);
138 
139  double ExtraInfo() const;
140 
141 protected:
142 
143  virtual std::string Name() const;
144 
145  virtual Singlet GetConstituentJet() const = 0;
146 
147  virtual boca::Jet GetJet() const = 0;
148 
149  virtual std::vector<boca::Jet> Jets() const = 0;
150 
151  virtual std::vector<LorentzVector<Momentum>> LorentzVectors() const = 0;
152 
154 
155 private:
156 
157  ClosestLepton closest_lepton_;
158 
159  Mutable<boca::Singlet> constituent_jet_;
160 
161  Mutable<boca::Jet> jet_;
162 
163  double extra_info_;
164 
165 };
166 
167 }
void SetExtraInfo(double extra_info)
Definition: Multiplet.cpp:102
Jet.
Definition: Jet.hh:15
Vector2< Angle > Angles(bool wrap_phi=false) const
Vector of rapidity and azimuth .
Definition: Multiplet.cpp:51
Vector2< Angle > AnglesMinTo(PseudoJet const &jet) const
Vector of rapidity and azimuth with minimal distance to jet.
Definition: PseudoJet.cpp:189
std::vector< boca::Jet > Constituents() const
All constituents.
Definition: Multiplet.cpp:66
boca::Mass Mass() const
Definition: Multiplet.cpp:61
boost::units::quantity< boost::units::si::plane_angle > Angle
Angle measured in radian.
Definition: Si.hh:35
boca::Singlet ConstituentJet() const
Jet of all constituents.
Definition: Multiplet.cpp:22
Lazy caching of variables.
Definition: Mutable.hh:19
Vector2< Angle > DeltaTo(Multiplet_ const &multiplet) const
Angular distance to a jet.
Definition: Multiplet.hh:88
boca::Jet Jet() const
Definition: Multiplet.cpp:29
Mutable< boca::EventShapes > event_shapes_
Definition: Multiplet.hh:153
Momentum Pt() const
Definition: Multiplet.cpp:36
Angle Phi() const
Azimuth .
Definition: Multiplet.cpp:46
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
Wrapper for a Jet in order to make it behave like a Multiplet.
Definition: Singlet.hh:19
double ExtraInfo() const
Definition: Multiplet.cpp:108
virtual boca::Jet GetJet() const =0
Angle DeltaRapTo(Multiplet_ const &multiplet) const
Difference to a jet.
Definition: Multiplet.hh:68
bool HasConstituents() const
Weather the jet as constituetns.
Definition: Multiplet.cpp:96
Calculate N-subjettiness of a jet.
Definition: SubJettiness.hh:17
virtual Singlet GetConstituentJet() const =0
Vector2< Angle > DeltaTo(PseudoJet const &jet) const
Angular distance to a jet.
Definition: PseudoJet.cpp:179
Angle DeltaPhiTo(Multiplet_ const &multiplet) const
Difference to a jet constrained to .
Definition: Multiplet.hh:58
boca::SubJettiness SubJettiness() const
Sub-jettiness.
Definition: Multiplet.cpp:56
Angle DeltaRapTo(PseudoJet const &jet) const
Difference to a jet.
Definition: PseudoJet.cpp:168
Angle DeltaRTo(Multiplet_ const &multiplet) const
Distance to a jet.
Definition: Multiplet.hh:78
ClosestLepton Lepton() const
Definition: Multiplet.cpp:17
Boosted Collider Analysis.
Definition: Analysis.hh:15
void SetClosestLepton(std::vector< boca::Lepton > const &leptons)
Definition: Multiplet.cpp:12
Multiplet base class
Definition: Multiplet.hh:21
virtual std::vector< LorentzVector< Momentum > > LorentzVectors() const =0
Angle Rap() const
Rapidity .
Definition: Multiplet.cpp:41
Angle DeltaPhiTo(PseudoJet const &jet) const
Difference to a jet constrained to .
Definition: PseudoJet.cpp:163
Definition: Identification.hh:10
Vector2< Angle > AnglesMinTo(Multiplet_ const &multiplet) const
Vector of rapidity and azimuth with minimal distance to jet.
Definition: Multiplet.hh:104
Definition: ClosestLepton.hh:10
Energy Momentum
Momentum measured in electronvolt.
Definition: ElectronVolt.hh:68
virtual std::vector< boca::Jet > Jets() const =0
Two dimensional Vector.
Definition: PseudoJet.hh:23
Energy Mass
Mass measured in electronvolt.
Definition: ElectronVolt.hh:62
virtual std::string Name() const
Definition: Multiplet.cpp:114