FCCAnalyses
Loading...
Searching...
No Matches
ReconstructedParticle.h
Go to the documentation of this file.
1#ifndef RECONSTRUCTEDPARTICLE_ANALYZERS_H
2#define RECONSTRUCTEDPARTICLE_ANALYZERS_H
3
4// Standard library
5#include <cmath>
6#include <vector>
7// ROOT
8#include "ROOT/RVec.hxx"
9#include "TLorentzVector.h"
10// EDM4hep
11#include "edm4hep/ReconstructedParticleData.h"
12#include "edm4hep/ParticleIDData.h"
13
14namespace FCCAnalyses ::ReconstructedParticle {
20 resonanceBuilder(float arg_resonance_mass);
22 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
23 operator()(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> legs);
24};
25
31 recoilBuilder(float arg_sqrts);
32 float m_sqrts = 240.0;
33 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator()(
34 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> inParticles);
35};
36
39 angular_separationBuilder( int arg_delta); // 0, 1, 2 = max, min, average
40 int m_delta = 0;
41 float operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in) ;
42 };
43
46 struct sel_type {
47 sel_type(const int type);
48 const int m_type;
49 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
50 operator()(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
51 };
52
55 struct sel_absType {
56 sel_absType(const int type);
57 const int m_type;
58 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
59 operator()(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
60 };
61
63 struct sel_pt {
64 sel_pt(float arg_min_pt);
65 float m_min_pt = 1.; //> transverse momentum threshold [GeV]
66 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
67 };
68
70 struct sel_eta {
71 sel_eta(float arg_min_eta);
72 float m_min_eta = 2.5; //> pseudorapidity threshold
73 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
74 };
75
77 struct sel_p {
78 sel_p(float arg_min_p, float arg_max_p = 1e10);
79 float m_min_p = 1.; //> momentum threshold [GeV]
80 float m_max_p = 1e10; //< momentum threshold [GeV]
81 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
82 };
83
85 struct sel_charge {
86 sel_charge(int arg_charge, bool arg_abs);
87 float m_charge; //> charge condition
88 bool m_abs;//> absolute value of the charge
89 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
90 };
91
93 struct sel_axis{
94 bool m_pos = 0; //> Which hemisphere to select, false/0=cosTheta<0 true/1=cosTheta>0
95 sel_axis(bool arg_pos);
96 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
97 operator()(ROOT::VecOps::RVec<float> angle,
98 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
99 };
100
102 struct sel_tag {
103 bool m_pass; // if pass is true, select tagged jets. Otherwise select anti-tagged ones
104 sel_tag(bool arg_pass);
105 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
106 operator()(ROOT::VecOps::RVec<bool> tags,
107 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
108 };
109
119 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
120 get(ROOT::VecOps::RVec<int> indexes,
121 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> inParticles);
122
124 ROOT::VecOps::RVec<float> get_pt(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
125
127 ROOT::VecOps::RVec<float> get_p(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
128
130 ROOT::VecOps::RVec<float> get_px(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
131
133 ROOT::VecOps::RVec<float> get_py(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
134
136 ROOT::VecOps::RVec<float> get_pz(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
137
139 ROOT::VecOps::RVec<float> get_eta(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
140
142 ROOT::VecOps::RVec<float> get_y(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
143
145 ROOT::VecOps::RVec<float> get_theta(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
146
148 ROOT::VecOps::RVec<float> get_phi(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
149
151 ROOT::VecOps::RVec<float> get_e(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
152
154 ROOT::VecOps::RVec<float> get_mass(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
155
157 ROOT::VecOps::RVec<float> get_charge(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
158
160 ROOT::VecOps::RVec<int> get_type(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
161
163 ROOT::VecOps::RVec<TLorentzVector> get_tlv(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
164
166 TLorentzVector get_tlv(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in, int index);
167
169 TLorentzVector get_tlv(edm4hep::ReconstructedParticleData in);
170
172 TLorentzVector get_P4vis(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
173
175 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
176 merge(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> x,
177 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> y);
178
180 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
181 remove(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> x,
182 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> y);
183
187 int get_n(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> inParticles);
188
192 ROOT::VecOps::RVec<bool>
193 getJet_btag(ROOT::VecOps::RVec<int> index,
194 ROOT::VecOps::RVec<edm4hep::ParticleIDData> pid,
195 ROOT::VecOps::RVec<float> values);
196
200 int getJet_ntags(ROOT::VecOps::RVec<bool> inBJetMask);
201 } // namespace FCCAnalyses::ReconstructedParticle
202
203#endif /* RECONSTRUCTEDPARTICLE_ANALYZERS_H */
TLorentzVector get_P4vis(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return visible 4-momentum vector
Definition ReconstructedParticle.cc:316
ROOT::VecOps::RVec< bool > getJet_btag(ROOT::VecOps::RVec< int > index, ROOT::VecOps::RVec< edm4hep::ParticleIDData > pid, ROOT::VecOps::RVec< float > values)
Returns the b-jet mask (vector of bools).
Definition ReconstructedParticle.cc:469
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > merge(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > x, ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > y)
concatenate both input vectors and return the resulting vector
Definition ReconstructedParticle.cc:265
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > remove(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > x, ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > y)
remove elements of vector y from vector x
Definition ReconstructedParticle.cc:275
int getJet_ntags(ROOT::VecOps::RVec< bool > inBJetMask)
Get number of b-jets.
Definition ReconstructedParticle.cc:482
FCC analyzers collection.
Definition Algorithms.h:15
return the angular separations (min / max / average) between a collection of particles
Definition ReconstructedParticle.h:38
Build the recoil from an arbitrary list of input ReconstructedPartilces and the center of mass energy...
Definition ReconstructedParticle.h:30
Build the two particle resonances from an arbitrary list of input ReconstructedPartilces.
Definition ReconstructedParticle.h:19
float m_resonance_mass
Definition ReconstructedParticle.h:21
select ReconstructedParticles by type absolute value Note: type might not correspond to PDG ID
Definition ReconstructedParticle.h:55
const int m_type
Definition ReconstructedParticle.h:57
select a list of reconstructed particles depending on the angle cosTheta axis
Definition ReconstructedParticle.h:93
select ReconstructedParticles with charge equal or in asolute value
Definition ReconstructedParticle.h:85
float m_charge
Definition ReconstructedParticle.h:87
bool m_abs
Definition ReconstructedParticle.h:88
select ReconstructedParticles with absolute pseudorapidity less than a maximum absolute value
Definition ReconstructedParticle.h:70
select ReconstructedParticles with momentum greater than a minimum value [GeV]
Definition ReconstructedParticle.h:77
select ReconstructedParticles with transverse momentum greater than a minimum value [GeV]
Definition ReconstructedParticle.h:63
select a list of reconstructed particles depending on the status of a certain boolean flag
Definition ReconstructedParticle.h:102
bool m_pass
Definition ReconstructedParticle.h:103
select ReconstructedParticles by type Note: type might not correspond to PDG ID
Definition ReconstructedParticle.h:46
const int m_type
Definition ReconstructedParticle.h:48