FCCAnalyses
Loading...
Searching...
No Matches
JetClusteringUtils.h
Go to the documentation of this file.
1
2#ifndef JETCLUSTERINGUTILS_ANALYZERS_H
3#define JETCLUSTERINGUTILS_ANALYZERS_H
4
5#include "Math/Vector4D.h"
6#include "ROOT/RVec.hxx"
7
9#include "fastjet/JetDefinition.hh"
10
11#include "TRandom3.h"
12
13namespace FCCAnalyses {
20 namespace JetClusteringUtils {
24 const int Nmax_dmerge = 10;
25
27 std::vector<fastjet::PseudoJet> set_pseudoJets(const ROOT::VecOps::RVec<float>& px,
28 const ROOT::VecOps::RVec<float>& py,
29 const ROOT::VecOps::RVec<float>& pz,
30 const ROOT::VecOps::RVec<float>& e);
31
40 std::vector<fastjet::PseudoJet> set_pseudoJets_xyzm(const ROOT::VecOps::RVec<float>& px,
41 const ROOT::VecOps::RVec<float>& py,
42 const ROOT::VecOps::RVec<float>& pz,
43 const ROOT::VecOps::RVec<float>& m);
44
46 ROOT::VecOps::RVec<fastjet::PseudoJet> get_pseudoJets(const JetClustering::FCCAnalysesJet& in);
47
49 std::vector<std::vector<int>> get_constituents(const JetClustering::FCCAnalysesJet& in);
50
53
55
57 ROOT::VecOps::RVec<float> get_px(const ROOT::VecOps::RVec<fastjet::PseudoJet>& in);
58
60 ROOT::VecOps::RVec<float> get_py(const ROOT::VecOps::RVec<fastjet::PseudoJet>& in);
61
63 ROOT::VecOps::RVec<float> get_pz(const ROOT::VecOps::RVec<fastjet::PseudoJet>& in);
64
66 ROOT::VecOps::RVec<float> get_e(const ROOT::VecOps::RVec<fastjet::PseudoJet>& in);
67
69 ROOT::VecOps::RVec<float> get_pt(const ROOT::VecOps::RVec<fastjet::PseudoJet>& in);
70
72 ROOT::VecOps::RVec<float> get_p(const ROOT::VecOps::RVec<fastjet::PseudoJet>& in);
73
75 ROOT::VecOps::RVec<float> get_m(const ROOT::VecOps::RVec<fastjet::PseudoJet>& in);
76
78 ROOT::VecOps::RVec<float> get_eta(const ROOT::VecOps::RVec<fastjet::PseudoJet>& in);
79
81 ROOT::VecOps::RVec<float> get_phi(const ROOT::VecOps::RVec<fastjet::PseudoJet>& in);
82
84 ROOT::VecOps::RVec<float> get_phi_std(const ROOT::VecOps::RVec<fastjet::PseudoJet>& in);
85
86
88 ROOT::VecOps::RVec<float> get_theta(const ROOT::VecOps::RVec<fastjet::PseudoJet>& in);
89
91 struct sel_pt {
92 sel_pt(float arg_min_pt);
93 float m_min_pt = 1.; //> transverse momentum threshold [GeV]
94 ROOT::VecOps::RVec<fastjet::PseudoJet> operator() (ROOT::VecOps::RVec<fastjet::PseudoJet> in);
95 };
96
99
100 JetClustering::FCCAnalysesJet build_FCCAnalysesJet(const std::vector<fastjet::PseudoJet>& in,
101 const std::vector<float>& dmerge,
102 const std::vector<float>& dmerge_max);
103
104 std::vector<fastjet::PseudoJet> build_jets(fastjet::ClusterSequence& cs, int exclusive, float cut, int sorted);
105
106 bool check(unsigned int n, int exclusive, float cut);
107
108 fastjet::RecombinationScheme recomb_scheme(int recombination);
109
110 std::vector<float> exclusive_dmerge(fastjet::ClusterSequence& cs, int do_dmarge_max);
111
112 // build the resonance from 2 <fastjet::PseudoJet> objects. Keep the closest to the mass given as input
115 resonanceBuilder(float arg_resonance_mass);
116 ROOT::VecOps::RVec<fastjet::PseudoJet> operator()(ROOT::VecOps::RVec<fastjet::PseudoJet> legs);
117 };
118
120 recoilBuilder(float arg_sqrts);
121 float m_sqrts = 240.0;
122 double operator() (ROOT::VecOps::RVec<fastjet::PseudoJet> in);
123 };
124 } // namespace JetClusteringUtils
125} // namespace FCCAnalyses
126
127#endif
ROOT::VecOps::RVec< float > get_phi(const ROOT::VecOps::RVec< fastjet::PseudoJet > &in)
Get jet phi.
Definition JetClusteringUtils.cc:135
ROOT::VecOps::RVec< float > get_m(const ROOT::VecOps::RVec< fastjet::PseudoJet > &in)
Get jet mass.
Definition JetClusteringUtils.cc:117
std::vector< fastjet::PseudoJet > set_pseudoJets_xyzm(const ROOT::VecOps::RVec< float > &px, const ROOT::VecOps::RVec< float > &py, const ROOT::VecOps::RVec< float > &pz, const ROOT::VecOps::RVec< float > &m)
Set fastjet pseudoJet for later reconstruction using px, py, pz and m.
Definition JetClusteringUtils.cc:44
std::vector< fastjet::PseudoJet > build_jets(fastjet::ClusterSequence &cs, int exclusive, float cut, int sorted)
Definition JetClusteringUtils.cc:215
ROOT::VecOps::RVec< fastjet::PseudoJet > get_pseudoJets(const JetClustering::FCCAnalysesJet &in)
Get fastjet pseudoJet after reconstruction from FCCAnalyses jets.
Definition JetClusteringUtils.cc:8
float get_exclusive_dmerge(const JetClustering::FCCAnalysesJet &in, int n)
return the dmin corresponding to the recombination that went from n+1 to n jets
Definition JetClusteringUtils.cc:17
JetClustering::FCCAnalysesJet build_FCCAnalysesJet(const std::vector< fastjet::PseudoJet > &in, const std::vector< float > &dmerge, const std::vector< float > &dmerge_max)
Definition JetClusteringUtils.cc:195
bool check(unsigned int n, int exclusive, float cut)
Definition JetClusteringUtils.cc:260
ROOT::VecOps::RVec< float > get_p(const ROOT::VecOps::RVec< fastjet::PseudoJet > &in)
Get jet p.
Definition JetClusteringUtils.cc:108
std::vector< float > exclusive_dmerge(fastjet::ClusterSequence &cs, int do_dmarge_max)
Definition JetClusteringUtils.cc:244
ROOT::VecOps::RVec< float > get_theta(const ROOT::VecOps::RVec< fastjet::PseudoJet > &in)
Get jet theta.
Definition JetClusteringUtils.cc:153
ROOT::VecOps::RVec< float > get_px(const ROOT::VecOps::RVec< fastjet::PseudoJet > &in)
Get jet px.
Definition JetClusteringUtils.cc:63
ROOT::VecOps::RVec< float > get_pz(const ROOT::VecOps::RVec< fastjet::PseudoJet > &in)
Get jet pz.
Definition JetClusteringUtils.cc:81
ROOT::VecOps::RVec< float > get_pt(const ROOT::VecOps::RVec< fastjet::PseudoJet > &in)
Get jet pt.
Definition JetClusteringUtils.cc:99
const int Nmax_dmerge
Maximum number of d_{n, n+1} that are kept in FCCAnalysesJet.
Definition JetClusteringUtils.h:24
ROOT::VecOps::RVec< float > get_phi_std(const ROOT::VecOps::RVec< fastjet::PseudoJet > &in)
Get jet phi.
Definition JetClusteringUtils.cc:144
JetClustering::FCCAnalysesJet initialise_FCCAnalysesJet()
Internal methods.
Definition JetClusteringUtils.cc:175
ROOT::VecOps::RVec< float > get_eta(const ROOT::VecOps::RVec< fastjet::PseudoJet > &in)
Get jet eta.
Definition JetClusteringUtils.cc:126
std::vector< std::vector< int > > get_constituents(const JetClustering::FCCAnalysesJet &in)
Get fastjet constituents after reconstruction from FCCAnalyses jets.
Definition JetClusteringUtils.cc:13
fastjet::RecombinationScheme recomb_scheme(int recombination)
Definition JetClusteringUtils.cc:266
ROOT::VecOps::RVec< float > get_py(const ROOT::VecOps::RVec< fastjet::PseudoJet > &in)
Get jet py.
Definition JetClusteringUtils.cc:72
ROOT::VecOps::RVec< float > get_e(const ROOT::VecOps::RVec< fastjet::PseudoJet > &in)
Get jet energy.
Definition JetClusteringUtils.cc:90
std::vector< fastjet::PseudoJet > set_pseudoJets(const ROOT::VecOps::RVec< float > &px, const ROOT::VecOps::RVec< float > &py, const ROOT::VecOps::RVec< float > &pz, const ROOT::VecOps::RVec< float > &e)
Set fastjet pseudoJet for later reconstruction.
Definition JetClusteringUtils.cc:31
float get_exclusive_dmerge_max(const JetClustering::FCCAnalysesJet &in, int n)
Definition JetClusteringUtils.cc:24
FCC analyzers collection.
Definition Algorithms.h:15
Definition JetClusteringUtils.h:119
recoilBuilder(float arg_sqrts)
Definition JetClusteringUtils.cc:327
float m_sqrts
Definition JetClusteringUtils.h:121
double operator()(ROOT::VecOps::RVec< fastjet::PseudoJet > in)
Definition JetClusteringUtils.cc:328
Definition JetClusteringUtils.h:113
float m_resonance_mass
Definition JetClusteringUtils.h:114
resonanceBuilder(float arg_resonance_mass)
Definition JetClusteringUtils.cc:289
ROOT::VecOps::RVec< fastjet::PseudoJet > operator()(ROOT::VecOps::RVec< fastjet::PseudoJet > legs)
Definition JetClusteringUtils.cc:293
Select clustered jets with transverse momentum greader than a minimum value [GeV].
Definition JetClusteringUtils.h:91
sel_pt(float arg_min_pt)
Definition JetClusteringUtils.cc:161
float m_min_pt
Definition JetClusteringUtils.h:93
ROOT::VecOps::RVec< fastjet::PseudoJet > operator()(ROOT::VecOps::RVec< fastjet::PseudoJet > in)
Definition JetClusteringUtils.cc:163
Structure to keep useful informations for the jets.
Definition JetClustering.h:22