FCCAnalyses
Loading...
Searching...
No Matches
ReconstructedParticle.h
Go to the documentation of this file.
1
2#ifndef RECONSTRUCTEDPARTICLE_ANALYZERS_H
3#define RECONSTRUCTEDPARTICLE_ANALYZERS_H
4
5// STL
6#include <cmath>
7#include <vector>
8
9// ROOT
10#include "TLorentzVector.h"
11#include "ROOT/RVec.hxx"
12
13// EDM4hep
14#include "edm4hep/ReconstructedParticleData.h"
15#include "edm4hep/ParticleIDData.h"
16
17namespace FCCAnalyses{
18
19namespace ReconstructedParticle{
20
24 resonanceBuilder(float arg_resonance_mass);
25 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator()(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> legs);
26 };
27
30 recoilBuilder(float arg_sqrts);
31 float m_sqrts = 240.0;
32 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in) ;
33 };
34
37 angular_separationBuilder( int arg_delta); // 0, 1, 2 = max, min, average
38 int m_delta = 0;
39 float operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in) ;
40 };
41
44 struct sel_type {
45 sel_type(const int type);
46 const int m_type;
47 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
48 operator()(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
49 };
50
53 struct sel_absType {
54 sel_absType(const int type);
55 const int m_type;
56 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
57 operator()(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
58 };
59
61 struct sel_pt {
62 sel_pt(float arg_min_pt);
63 float m_min_pt = 1.; //> transverse momentum threshold [GeV]
64 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
65 };
66
68 struct sel_eta {
69 sel_eta(float arg_min_eta);
70 float m_min_eta = 2.5; //> pseudorapidity threshold
71 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
72 };
73
75 struct sel_p {
76 sel_p(float arg_min_p, float arg_max_p = 1e10);
77 float m_min_p = 1.; //> momentum threshold [GeV]
78 float m_max_p = 1e10; //< momentum threshold [GeV]
79 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
80 };
81
83 struct sel_charge {
84 sel_charge(int arg_charge, bool arg_abs);
85 float m_charge; //> charge condition
86 bool m_abs;//> absolute value of the charge
87 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
88 };
89
91 struct sel_axis{
92 bool m_pos = 0; //> Which hemisphere to select, false/0=cosTheta<0 true/1=cosTheta>0
93 sel_axis(bool arg_pos);
94 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator()(ROOT::VecOps::RVec<float> angle, ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
95 };
96
98 struct sel_tag {
99 bool m_pass; // if pass is true, select tagged jets. Otherwise select anti-tagged ones
100 sel_tag(bool arg_pass);
101 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator() (ROOT::VecOps::RVec<bool> tags, ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
102 };
103
105 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> get(ROOT::VecOps::RVec<int> index, ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
106
108 ROOT::VecOps::RVec<float> get_pt(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
109
111 ROOT::VecOps::RVec<float> get_p(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
112
114 ROOT::VecOps::RVec<float> get_px(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
115
117 ROOT::VecOps::RVec<float> get_py(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
118
120 ROOT::VecOps::RVec<float> get_pz(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
121
123 ROOT::VecOps::RVec<float> get_eta(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
124
126 ROOT::VecOps::RVec<float> get_y(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
127
129 ROOT::VecOps::RVec<float> get_theta(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
130
132 ROOT::VecOps::RVec<float> get_phi(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
133
135 ROOT::VecOps::RVec<float> get_e(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
136
138 ROOT::VecOps::RVec<float> get_mass(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
139
141 ROOT::VecOps::RVec<float> get_charge(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
142
144 ROOT::VecOps::RVec<int> get_type(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
145
147 ROOT::VecOps::RVec<TLorentzVector> get_tlv(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
148
150 TLorentzVector get_tlv(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in, int index);
151
153 TLorentzVector get_tlv(edm4hep::ReconstructedParticleData in);
154
156 TLorentzVector get_P4vis(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
157
159 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> merge(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> x, ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> y);
160
162 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> remove( ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> x, ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> y);
163
165 int get_n(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
166
168 ROOT::VecOps::RVec<bool> getJet_btag(ROOT::VecOps::RVec<int> index, ROOT::VecOps::RVec<edm4hep::ParticleIDData> pid, ROOT::VecOps::RVec<float> values);
169
171 int getJet_ntags(ROOT::VecOps::RVec<bool> in);
172}//end NS ReconstructedParticle
173
174}//end NS FCCAnalyses
175#endif
TLorentzVector get_P4vis(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return visible 4-momentum vector
Definition ReconstructedParticle.cc:306
ROOT::VecOps::RVec< float > get_pz(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the momenta of the input ReconstructedParticles
Definition ReconstructedParticle.cc:380
ROOT::VecOps::RVec< float > get_px(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the momenta of the input ReconstructedParticles
Definition ReconstructedParticle.cc:363
ROOT::VecOps::RVec< float > get_y(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the rapidity of the input ReconstructedParticles
Definition ReconstructedParticle.cc:396
ROOT::VecOps::RVec< float > get_phi(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the phi of the input ReconstructedParticles
Definition ReconstructedParticle.cc:335
ROOT::VecOps::RVec< float > get_mass(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the masses of the input ReconstructedParticles
Definition ReconstructedParticle.cc:317
ROOT::VecOps::RVec< float > get_p(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the momenta of the input ReconstructedParticles
Definition ReconstructedParticle.cc:353
ROOT::VecOps::RVec< TLorentzVector > get_tlv(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the TlorentzVector of the input ReconstructedParticles
Definition ReconstructedParticle.cc:416
ROOT::VecOps::RVec< float > get_eta(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the pseudo-rapidity of the input ReconstructedParticles
Definition ReconstructedParticle.cc:325
ROOT::VecOps::RVec< float > get_py(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the momenta of the input ReconstructedParticles
Definition ReconstructedParticle.cc:372
ROOT::VecOps::RVec< float > get_pt(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the transverse momenta of the input ReconstructedParticles
Definition ReconstructedParticle.cc:243
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:251
ROOT::VecOps::RVec< float > get_theta(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the theta of the input ReconstructedParticles
Definition ReconstructedParticle.cc:406
ROOT::VecOps::RVec< float > get_e(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the energy of the input ReconstructedParticles
Definition ReconstructedParticle.cc:345
ROOT::VecOps::RVec< bool > getJet_btag(ROOT::VecOps::RVec< int > index, ROOT::VecOps::RVec< edm4hep::ParticleIDData > pid, ROOT::VecOps::RVec< float > values)
returns the bjet flavour
Definition ReconstructedParticle.cc:459
int get_n(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the size of the input collection
Definition ReconstructedParticle.cc:453
int getJet_ntags(ROOT::VecOps::RVec< bool > in)
get number of b-jets
Definition ReconstructedParticle.cc:473
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > get(ROOT::VecOps::RVec< int > index, ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return reconstructed particles
Definition ReconstructedParticle.cc:295
ROOT::VecOps::RVec< int > get_type(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the type of the input ReconstructedParticles
Definition ReconstructedParticle.cc:440
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:261
ROOT::VecOps::RVec< float > get_charge(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
return the charges of the input ReconstructedParticles
Definition ReconstructedParticle.cc:388
FCC analyzers collection.
Definition Algorithms.h:15
return the angular separations (min / max / average) between a collection of particles
Definition ReconstructedParticle.h:36
angular_separationBuilder(int arg_delta)
Definition ReconstructedParticle.cc:209
float operator()(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
Definition ReconstructedParticle.cc:210
int m_delta
Definition ReconstructedParticle.h:38
build the recoil from an arbitrary list of input ReconstructedPartilces and the center of mass energy
Definition ReconstructedParticle.h:29
float m_sqrts
Definition ReconstructedParticle.h:31
recoilBuilder(float arg_sqrts)
Definition ReconstructedParticle.cc:160
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > operator()(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
Definition ReconstructedParticle.cc:161
build the resonance from 2 particles from an arbitrary list of input ReconstructedPartilces....
Definition ReconstructedParticle.h:22
resonanceBuilder(float arg_resonance_mass)
Definition ReconstructedParticle.cc:122
float m_resonance_mass
Definition ReconstructedParticle.h:23
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > operator()(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > legs)
Definition ReconstructedParticle.cc:123
select ReconstructedParticles by type absolute value Note: type might not correspond to PDG ID
Definition ReconstructedParticle.h:53
const int m_type
Definition ReconstructedParticle.h:55
sel_absType(const int type)
sel_absType
Definition ReconstructedParticle.cc:38
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > operator()(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
Definition ReconstructedParticle.cc:45
select a list of reconstructed particles depending on the angle cosTheta axis
Definition ReconstructedParticle.h:91
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > operator()(ROOT::VecOps::RVec< float > angle, ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
Definition ReconstructedParticle.cc:180
sel_axis(bool arg_pos)
Definition ReconstructedParticle.cc:179
bool m_pos
Definition ReconstructedParticle.h:92
select ReconstructedParticles with charge equal or in asolute value
Definition ReconstructedParticle.h:83
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > operator()(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
Definition ReconstructedParticle.cc:110
float m_charge
Definition ReconstructedParticle.h:85
bool m_abs
Definition ReconstructedParticle.h:86
sel_charge(int arg_charge, bool arg_abs)
Definition ReconstructedParticle.cc:108
select ReconstructedParticles with absolute pseudorapidity less than a maximum absolute value
Definition ReconstructedParticle.h:68
sel_eta(float arg_min_eta)
Definition ReconstructedParticle.cc:76
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > operator()(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
Definition ReconstructedParticle.cc:77
float m_min_eta
Definition ReconstructedParticle.h:70
select ReconstructedParticles with momentum greater than a minimum value [GeV]
Definition ReconstructedParticle.h:75
float m_min_p
Definition ReconstructedParticle.h:77
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > operator()(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
Definition ReconstructedParticle.cc:93
sel_p(float arg_min_p, float arg_max_p=1e10)
Definition ReconstructedParticle.cc:92
float m_max_p
Definition ReconstructedParticle.h:78
select ReconstructedParticles with transverse momentum greater than a minimum value [GeV]
Definition ReconstructedParticle.h:61
float m_min_pt
Definition ReconstructedParticle.h:63
sel_pt(float arg_min_pt)
sel_pt
Definition ReconstructedParticle.cc:63
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > operator()(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
Definition ReconstructedParticle.cc:64
select a list of reconstructed particles depending on the status of a certain boolean flag
Definition ReconstructedParticle.h:98
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > operator()(ROOT::VecOps::RVec< bool > tags, ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
Definition ReconstructedParticle.cc:191
sel_tag(bool arg_pass)
Definition ReconstructedParticle.cc:190
bool m_pass
Definition ReconstructedParticle.h:99
select ReconstructedParticles by type Note: type might not correspond to PDG ID
Definition ReconstructedParticle.h:44
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > operator()(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > in)
Definition ReconstructedParticle.cc:21
const int m_type
Definition ReconstructedParticle.h:46
sel_type(const int type)
sel_type
Definition ReconstructedParticle.cc:19