1#ifndef RECONSTRUCTEDPARTICLE_ANALYZERS_H
2#define RECONSTRUCTEDPARTICLE_ANALYZERS_H
8#include "ROOT/RVec.hxx"
9#include "TLorentzVector.h"
11#include "edm4hep/ReconstructedParticleData.h"
12#include "edm4hep/ParticleIDData.h"
22 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
23 operator()(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> legs);
32 float m_sqrts = 240.0;
33 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator()(
34 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> inParticles);
41 float operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in) ;
49 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
50 operator()(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
58 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
59 operator()(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
66 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
72 float m_min_eta = 2.5;
73 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
78 sel_p(
float arg_min_p,
float arg_max_p = 1e10);
81 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
89 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> operator() (ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
96 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
97 operator()(ROOT::VecOps::RVec<float> angle,
98 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
105 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
106 operator()(ROOT::VecOps::RVec<bool> tags,
107 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
119 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
120 get(ROOT::VecOps::RVec<int> indexes,
121 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> inParticles);
124 ROOT::VecOps::RVec<float> get_pt(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
127 ROOT::VecOps::RVec<float> get_p(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
130 ROOT::VecOps::RVec<float> get_px(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
133 ROOT::VecOps::RVec<float> get_py(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
136 ROOT::VecOps::RVec<float> get_pz(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
139 ROOT::VecOps::RVec<float> get_eta(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
142 ROOT::VecOps::RVec<float> get_y(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
145 ROOT::VecOps::RVec<float> get_theta(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
148 ROOT::VecOps::RVec<float> get_phi(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
151 ROOT::VecOps::RVec<float> get_e(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
154 ROOT::VecOps::RVec<float> get_mass(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
157 ROOT::VecOps::RVec<float> get_charge(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
160 ROOT::VecOps::RVec<int> get_type(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
163 ROOT::VecOps::RVec<TLorentzVector> get_tlv(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
166 TLorentzVector get_tlv(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in,
int index);
169 TLorentzVector get_tlv(edm4hep::ReconstructedParticleData in);
172 TLorentzVector
get_P4vis(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> in);
175 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
176 merge(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> x,
177 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> y);
180 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
181 remove(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> x,
182 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> y);
187 int get_n(ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> inParticles);
192 ROOT::VecOps::RVec<bool>
194 ROOT::VecOps::RVec<edm4hep::ParticleIDData> pid,
195 ROOT::VecOps::RVec<float> values);
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