FCCAnalyses
Loading...
Searching...
No Matches
JetConstituentsUtils.h
Go to the documentation of this file.
1#ifndef FCCAnalyses_JetConstituentsUtils_h
2#define FCCAnalyses_JetConstituentsUtils_h
3
4// ROOT
5#include "ROOT/RVec.hxx"
6#include "TLorentzVector.h"
7#include "TMath.h"
8#include "TRotation.h"
9#include "TVector3.h"
10// EDM4hep
11#include "edm4hep/CalorimeterHitData.h"
12#include "edm4hep/ClusterData.h"
13#include "edm4hep/MCParticleData.h"
14#include "edm4hep/RecDqdxData.h"
15#include "edm4hep/ReconstructedParticleData.h"
16#if __has_include("edm4hep/TrackerHit3DData.h")
17#include "edm4hep/TrackerHit3DData.h"
18#else
19#include "edm4hep/TrackerHitData.h"
20namespace edm4hep {
21 using TrackerHit3DData = edm4hep::TrackerHitData;
22}
23#endif
24#include "edm4hep/TrackData.h"
25#include "edm4hep/TrackState.h"
26// FastJet
27#include "fastjet/JetDefinition.hh"
28// FCCAnalyses
30
31namespace FCCAnalyses {
32 namespace JetConstituentsUtils {
33 namespace rv = ROOT::VecOps;
34 using FCCAnalysesJetConstituents = rv::RVec<edm4hep::ReconstructedParticleData>;
35 using FCCAnalysesJetConstituentsData = rv::RVec<float>;
36
38 rv::RVec<FCCAnalysesJetConstituents> build_constituents(const rv::RVec<edm4hep::ReconstructedParticleData>&,
39 const rv::RVec<edm4hep::ReconstructedParticleData>&);
40
41 rv::RVec<FCCAnalysesJetConstituents> build_constituents_cluster(const rv::RVec<edm4hep::ReconstructedParticleData>& rps,
42 const std::vector<std::vector<int>>& indices);
43
45 FCCAnalysesJetConstituents get_jet_constituents(const rv::RVec<FCCAnalysesJetConstituents>&, int);
47 rv::RVec<FCCAnalysesJetConstituents> get_constituents(const rv::RVec<FCCAnalysesJetConstituents>&,
48 const rv::RVec<int>&);
49
50
51 //sorting jets
52 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> jets_sorting_on_nconst(const rv::RVec<edm4hep::ReconstructedParticleData>&);
53 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> jets_sorting_on_energy(const rv::RVec<edm4hep::ReconstructedParticleData>&);
54
55
56 rv::RVec<FCCAnalysesJetConstituentsData> get_Bz(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
57 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
58
59 rv::RVec<FCCAnalysesJetConstituentsData> get_pt(const rv::RVec<FCCAnalysesJetConstituents>&);
60 rv::RVec<FCCAnalysesJetConstituentsData> get_p(const rv::RVec<FCCAnalysesJetConstituents>&);
61 rv::RVec<FCCAnalysesJetConstituentsData> get_e(const rv::RVec<FCCAnalysesJetConstituents>&);
62 rv::RVec<FCCAnalysesJetConstituentsData> get_theta(const rv::RVec<FCCAnalysesJetConstituents>&);
63 rv::RVec<FCCAnalysesJetConstituentsData> get_phi(const rv::RVec<FCCAnalysesJetConstituents>&);
64 rv::RVec<FCCAnalysesJetConstituentsData> get_type(const rv::RVec<FCCAnalysesJetConstituents>&);
65 rv::RVec<FCCAnalysesJetConstituentsData> get_charge(const rv::RVec<FCCAnalysesJetConstituents>&);
66
67 //displacement
68 rv::RVec<FCCAnalysesJetConstituentsData> get_d0(const rv::RVec<FCCAnalysesJetConstituents>&,
69 const ROOT::VecOps::RVec<edm4hep::TrackState>&);
70
71 rv::RVec<FCCAnalysesJetConstituentsData> get_z0(const rv::RVec<FCCAnalysesJetConstituents>& ,
72 const ROOT::VecOps::RVec<edm4hep::TrackState>&);
73
74 rv::RVec<FCCAnalysesJetConstituentsData> get_phi0(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
75 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
76
77 rv::RVec<FCCAnalysesJetConstituentsData> get_omega(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
78 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
79
80 rv::RVec<FCCAnalysesJetConstituentsData> get_tanLambda(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
81 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
82
83
84 rv::RVec<FCCAnalysesJetConstituentsData> XPtoPar_dxy(const rv::RVec<FCCAnalysesJetConstituents>&,
85 const ROOT::VecOps::RVec<edm4hep::TrackState>&,
86 const TLorentzVector& V, // primary vertex
87 const float&);
88 rv::RVec<FCCAnalysesJetConstituentsData> XPtoPar_dz(const rv::RVec<FCCAnalysesJetConstituents>&,
89 const ROOT::VecOps::RVec<edm4hep::TrackState>&,
90 const TLorentzVector& V, // primary vertex
91 const float&);
92 rv::RVec<FCCAnalysesJetConstituentsData> XPtoPar_phi(const rv::RVec<FCCAnalysesJetConstituents>&,
93 const ROOT::VecOps::RVec<edm4hep::TrackState>&,
94 const TLorentzVector& V, // primary vertex
95 const float&);
96 rv::RVec<FCCAnalysesJetConstituentsData> XPtoPar_C(const rv::RVec<FCCAnalysesJetConstituents>&,
97 const ROOT::VecOps::RVec<edm4hep::TrackState>&,
98 const float&);
99 rv::RVec<FCCAnalysesJetConstituentsData> XPtoPar_ct(const rv::RVec<FCCAnalysesJetConstituents>&,
100 const ROOT::VecOps::RVec<edm4hep::TrackState>&,
101 const float&);
102
103 //covariance matrix
104 //diagonal
105 rv::RVec<FCCAnalysesJetConstituentsData> get_omega_cov(const rv::RVec<FCCAnalysesJetConstituents>&,
106 const ROOT::VecOps::RVec<edm4hep::TrackState>&);
107
108 rv::RVec<FCCAnalysesJetConstituentsData> get_d0_cov(const rv::RVec<FCCAnalysesJetConstituents>&,
109 const ROOT::VecOps::RVec<edm4hep::TrackState>& );
110
111 rv::RVec<FCCAnalysesJetConstituentsData> get_z0_cov(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
112 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
113
114 rv::RVec<FCCAnalysesJetConstituentsData> get_phi0_cov(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
115 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
116
117 rv::RVec<FCCAnalysesJetConstituentsData> get_tanlambda_cov(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
118 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
119 //off-diag
120 rv::RVec<FCCAnalysesJetConstituentsData> get_d0_z0_cov(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
121 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
122
123 rv::RVec<FCCAnalysesJetConstituentsData> get_phi0_d0_cov(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
124 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
125
126 rv::RVec<FCCAnalysesJetConstituentsData> get_phi0_z0_cov(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
127 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
128
129 rv::RVec<FCCAnalysesJetConstituentsData> get_tanlambda_phi0_cov(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
130 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
131
132 rv::RVec<FCCAnalysesJetConstituentsData> get_tanlambda_d0_cov(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
133 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
134
135 rv::RVec<FCCAnalysesJetConstituentsData> get_tanlambda_z0_cov(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
136 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
137
138 rv::RVec<FCCAnalysesJetConstituentsData> get_omega_tanlambda_cov(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
139 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
140
141 rv::RVec<FCCAnalysesJetConstituentsData> get_omega_phi0_cov(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
142 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
143
144 rv::RVec<FCCAnalysesJetConstituentsData> get_omega_d0_cov(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
145 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
146
147 rv::RVec<FCCAnalysesJetConstituentsData> get_omega_z0_cov(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
148 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
149
165 rv::RVec<FCCAnalysesJetConstituentsData> get_dndx(
166 const rv::RVec<FCCAnalysesJetConstituents> &jetConstituents,
167 const TrackUtils::TrackDqdxHandler &dNdxHandler,
168 const rv::RVec<edm4hep::TrackData> &trackColl,
169 const rv::RVec<FCCAnalysesJetConstituentsData> isJetConstChargedHad);
170
171 rv::RVec<FCCAnalysesJetConstituentsData> get_Sip2dVal(const rv::RVec<edm4hep::ReconstructedParticleData>& jets,
172 const rv::RVec<FCCAnalysesJetConstituents>& jcs,
173 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
174
175 rv::RVec<FCCAnalysesJetConstituentsData> get_Sip2dVal_cluster(const rv::RVec<fastjet::PseudoJet>& jets,
176 const rv::RVec<FCCAnalysesJetConstituents>& jcs,
177 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
178
179
180 rv::RVec<FCCAnalysesJetConstituentsData> get_Sip2dVal_clusterV(const rv::RVec<fastjet::PseudoJet>& jets,
181 const rv::RVec<FCCAnalysesJetConstituentsData>& D0,
182 const rv::RVec<FCCAnalysesJetConstituentsData>& phi0,
183 const float Bz);
184
185
186 rv::RVec<FCCAnalysesJetConstituentsData> get_Sip2dSig(const rv::RVec<FCCAnalysesJetConstituentsData>& Sip2dVals,
187 const rv::RVec<FCCAnalysesJetConstituentsData>& err2_D0);
188
189 rv::RVec<FCCAnalysesJetConstituentsData> get_Sip3dVal(const rv::RVec<edm4hep::ReconstructedParticleData>& jets,
190 const rv::RVec<FCCAnalysesJetConstituents>& jcs,
191 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
192
193
194 rv::RVec<FCCAnalysesJetConstituentsData> get_Sip3dVal_cluster(const rv::RVec<fastjet::PseudoJet>& jets,
195 const rv::RVec<FCCAnalysesJetConstituents>& jcs,
196 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
197
198 rv::RVec<FCCAnalysesJetConstituentsData> get_Sip3dVal_clusterV(const rv::RVec<fastjet::PseudoJet>& jets,
199 const rv::RVec<FCCAnalysesJetConstituentsData>& D0,
200 const rv::RVec<FCCAnalysesJetConstituentsData>& Z0,
201 const rv::RVec<FCCAnalysesJetConstituentsData>& phi0,
202 const float Bz);
203
204 rv::RVec<FCCAnalysesJetConstituentsData> get_Sip3dSig(const rv::RVec<FCCAnalysesJetConstituentsData>& Sip3dVals,
205 const rv::RVec<FCCAnalysesJetConstituentsData>& err2_D0,
206 const rv::RVec<FCCAnalysesJetConstituentsData>& err2_Z0);
207
208 rv::RVec<FCCAnalysesJetConstituentsData> get_JetDistVal(const rv::RVec<edm4hep::ReconstructedParticleData>& jets,
209 const rv::RVec<FCCAnalysesJetConstituents>& jcs,
210 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
211
212 rv::RVec<FCCAnalysesJetConstituentsData> get_JetDistVal_cluster(const rv::RVec<fastjet::PseudoJet>& jets,
213 const rv::RVec<FCCAnalysesJetConstituents>& jcs,
214 const ROOT::VecOps::RVec<edm4hep::TrackState>& tracks);
215
216 rv::RVec<FCCAnalysesJetConstituentsData> get_JetDistVal_clusterV(const rv::RVec<fastjet::PseudoJet>& jets,
217 const rv::RVec<FCCAnalysesJetConstituents>& jcs,
218 const rv::RVec<FCCAnalysesJetConstituentsData>& D0,
219 const rv::RVec<FCCAnalysesJetConstituentsData>& Z0,
220 const rv::RVec<FCCAnalysesJetConstituentsData>& phi0,
221 const float Bz);
222
223 rv::RVec<FCCAnalysesJetConstituentsData> get_JetDistSig(const rv::RVec<FCCAnalysesJetConstituentsData>& JetDistVal,
224 const rv::RVec<FCCAnalysesJetConstituentsData>& err2_D0,
225 const rv::RVec<FCCAnalysesJetConstituentsData>& err2_Z0);
226
227 rv::RVec<FCCAnalysesJetConstituentsData> get_mtof(const rv::RVec<FCCAnalysesJetConstituents>& jcs,
228 const rv::RVec<float>& track_L,
229 const rv::RVec<edm4hep::TrackData>& trackdata,
230 const rv::RVec<edm4hep::TrackerHit3DData>& trackerhits,
231 const rv::RVec<edm4hep::ClusterData>& gammadata,
232 const rv::RVec<edm4hep::ClusterData>& nhdata,
233 const rv::RVec<edm4hep::CalorimeterHitData>& calohits,
234 const TLorentzVector& V // primary vertex
235 );
236
237
238 rv::RVec<FCCAnalysesJetConstituentsData> get_PIDs(const ROOT::VecOps::RVec< int > recin,
239 const ROOT::VecOps::RVec< int > mcin,
240 const rv::RVec<edm4hep::ReconstructedParticleData>& RecPart,
241 const rv::RVec<edm4hep::MCParticleData>& Particle,
242 const rv::RVec<edm4hep::ReconstructedParticleData>& Jets);
243
244 rv::RVec<FCCAnalysesJetConstituentsData> get_PIDs_cluster(const ROOT::VecOps::RVec< int > recin,
245 const ROOT::VecOps::RVec< int > mcin,
246 const rv::RVec<edm4hep::ReconstructedParticleData>& RecPart,
247 const rv::RVec<edm4hep::MCParticleData>& Particle,
248 const std::vector<std::vector<int>>& indices);
249
250 rv::RVec<FCCAnalysesJetConstituentsData> get_isMu(const rv::RVec<FCCAnalysesJetConstituents>& jcs);
251 rv::RVec<FCCAnalysesJetConstituentsData> get_isEl(const rv::RVec<FCCAnalysesJetConstituents>& jcs);
252 rv::RVec<FCCAnalysesJetConstituentsData> get_isChargedHad(const rv::RVec<FCCAnalysesJetConstituents>& jcs);
253 rv::RVec<FCCAnalysesJetConstituentsData> get_isGamma(const rv::RVec<FCCAnalysesJetConstituents>& jcs);
254 rv::RVec<FCCAnalysesJetConstituentsData> get_isNeutralHad(const rv::RVec<FCCAnalysesJetConstituents>& jcs);
255
256 //countings
257 int count_jets(rv::RVec<FCCAnalysesJetConstituents> jets);
258 rv::RVec<int> count_consts(rv::RVec<FCCAnalysesJetConstituents> jets);
259 rv::RVec<int> count_type(const rv::RVec<FCCAnalysesJetConstituentsData>& isType);
260
261
262
263 rv::RVec<FCCAnalysesJetConstituentsData> get_erel(const rv::RVec<edm4hep::ReconstructedParticleData>& jets,
264 const rv::RVec<FCCAnalysesJetConstituents>& jcs);
265 rv::RVec<FCCAnalysesJetConstituentsData> get_erel_cluster(const rv::RVec<fastjet::PseudoJet>& jets,
266 const rv::RVec<FCCAnalysesJetConstituents>& jcs);
267
268 rv::RVec<FCCAnalysesJetConstituentsData> get_erel_log(const rv::RVec<edm4hep::ReconstructedParticleData>& jets,
269 const rv::RVec<FCCAnalysesJetConstituents>& jcs);
270 rv::RVec<FCCAnalysesJetConstituentsData> get_erel_log_cluster(const rv::RVec<fastjet::PseudoJet>& jets,
271 const rv::RVec<FCCAnalysesJetConstituents>& jcs);
272
273 rv::RVec<FCCAnalysesJetConstituentsData> get_thetarel(const rv::RVec<edm4hep::ReconstructedParticleData>& jets,
274 const rv::RVec<FCCAnalysesJetConstituents>& jcs);
275 rv::RVec<FCCAnalysesJetConstituentsData> get_thetarel_cluster(const rv::RVec<fastjet::PseudoJet>& jets,
276 const rv::RVec<FCCAnalysesJetConstituents>& jcs);
277
278 rv::RVec<FCCAnalysesJetConstituentsData> get_phirel(const rv::RVec<edm4hep::ReconstructedParticleData>& jets,
279 const rv::RVec<FCCAnalysesJetConstituents>& jcs);
280 rv::RVec<FCCAnalysesJetConstituentsData> get_phirel_cluster(const rv::RVec<fastjet::PseudoJet>& jets,
281 const rv::RVec<FCCAnalysesJetConstituents>& jcs);
282
283 //residues
284 rv::RVec<TLorentzVector> compute_tlv_jets(const rv::RVec<fastjet::PseudoJet>& jets);
285 rv::RVec<TLorentzVector> sum_tlv_constituents(const rv::RVec<FCCAnalysesJetConstituents>& jets);
286 float InvariantMass(const TLorentzVector& tlv1, const TLorentzVector& tlv2);
287
291 rv::RVec<double> all_invariant_masses(rv::RVec<TLorentzVector> AllJets);
292 rv::RVec<double> compute_residue_energy(const rv::RVec<TLorentzVector>& tlv_jet,
293 const rv::RVec<TLorentzVector>& sum_tlv_jcs);
294 rv::RVec<double> compute_residue_pt(const rv::RVec<TLorentzVector>& tlv_jet,
295 const rv::RVec<TLorentzVector>& sum_tlv_jcs);
296 rv::RVec<double> compute_residue_phi(const rv::RVec<TLorentzVector>& tlv_jet,
297 const rv::RVec<TLorentzVector>& sum_tlv_jcs);
298 rv::RVec<double> compute_residue_theta(const rv::RVec<TLorentzVector>& tlv_jet,
299 const rv::RVec<TLorentzVector>& sum_tlv_jcs);
300 rv::RVec<double> compute_residue_px(const rv::RVec<TLorentzVector>& tlv_jet, const rv::RVec<TLorentzVector>& sum_tlv_jcs);
301 rv::RVec<double> compute_residue_py(const rv::RVec<TLorentzVector>& tlv_jet, const rv::RVec<TLorentzVector>& sum_tlv_jcs);
302 rv::RVec<double> compute_residue_pz(const rv::RVec<TLorentzVector>& tlv_jet, const rv::RVec<TLorentzVector>& sum_tlv_jcs);
303
304 } // namespace JetConstituentsUtils
305} // namespace FCCAnalyses
306
307#endif
Adjusted utility class to invert the relations between RecDqdx and Track.
Definition TrackUtils.h:26
rv::RVec< FCCAnalysesJetConstituentsData > get_isMu(const rv::RVec< FCCAnalysesJetConstituents > &jcs)
Definition JetConstituentsUtils.cc:1112
rv::RVec< FCCAnalysesJetConstituentsData > get_isGamma(const rv::RVec< FCCAnalysesJetConstituents > &jcs)
Definition JetConstituentsUtils.cc:1185
rv::RVec< FCCAnalysesJetConstituentsData > get_erel_cluster(const rv::RVec< fastjet::PseudoJet > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs)
Definition JetConstituentsUtils.cc:920
rv::RVec< FCCAnalysesJetConstituentsData > get_phirel(const rv::RVec< edm4hep::ReconstructedParticleData > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs)
Definition JetConstituentsUtils.cc:991
rv::RVec< FCCAnalysesJetConstituentsData > get_thetarel(const rv::RVec< edm4hep::ReconstructedParticleData > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs)
Definition JetConstituentsUtils.cc:939
rv::RVec< FCCAnalysesJetConstituentsData > get_isNeutralHad(const rv::RVec< FCCAnalysesJetConstituents > &jcs)
Definition JetConstituentsUtils.cc:1160
rv::RVec< FCCAnalysesJetConstituentsData > get_z0(const rv::RVec< FCCAnalysesJetConstituents > &, const ROOT::VecOps::RVec< edm4hep::TrackState > &)
Definition JetConstituentsUtils.cc:201
rv::RVec< FCCAnalysesJetConstituentsData > get_JetDistVal(const rv::RVec< edm4hep::ReconstructedParticleData > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:621
rv::RVec< FCCAnalysesJetConstituentsData > get_theta(const rv::RVec< FCCAnalysesJetConstituents > &)
Definition JetConstituentsUtils.cc:141
rv::RVec< FCCAnalysesJetConstituentsData > get_pt(const rv::RVec< FCCAnalysesJetConstituents > &)
Definition JetConstituentsUtils.cc:126
float InvariantMass(const TLorentzVector &tlv1, const TLorentzVector &tlv2)
Definition JetConstituentsUtils.cc:1274
rv::RVec< FCCAnalysesJetConstituentsData > get_phi0_d0_cov(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:306
rv::RVec< TLorentzVector > sum_tlv_constituents(const rv::RVec< FCCAnalysesJetConstituents > &jets)
Definition JetConstituentsUtils.cc:1256
rv::RVec< FCCAnalysesJetConstituentsData > get_JetDistVal_clusterV(const rv::RVec< fastjet::PseudoJet > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs, const rv::RVec< FCCAnalysesJetConstituentsData > &D0, const rv::RVec< FCCAnalysesJetConstituentsData > &Z0, const rv::RVec< FCCAnalysesJetConstituentsData > &phi0, const float Bz)
Definition JetConstituentsUtils.cc:687
rv::RVec< double > all_invariant_masses(rv::RVec< TLorentzVector > AllJets)
all_invariant_masses takes an RVec of TLorentzVectors of jets and computes the invariant masses of al...
Definition JetConstituentsUtils.cc:1284
rv::RVec< TLorentzVector > compute_tlv_jets(const rv::RVec< fastjet::PseudoJet > &jets)
Definition JetConstituentsUtils.cc:1244
rv::RVec< FCCAnalysesJetConstituentsData > get_omega(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:213
rv::RVec< FCCAnalysesJetConstituentsData > get_phi0_cov(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:288
rv::RVec< edm4hep::ReconstructedParticleData > FCCAnalysesJetConstituents
Definition JetConstituentsUtils.h:34
rv::RVec< FCCAnalysesJetConstituentsData > get_tanlambda_z0_cov(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:330
rv::RVec< FCCAnalysesJetConstituents > build_constituents(const rv::RVec< edm4hep::ReconstructedParticleData > &, const rv::RVec< edm4hep::ReconstructedParticleData > &)
Build the collection of constituents (mapping jet -> reconstructed particles) for all jets in event.
Definition JetConstituentsUtils.cc:27
rv::RVec< FCCAnalysesJetConstituents > build_constituents_cluster(const rv::RVec< edm4hep::ReconstructedParticleData > &rps, const std::vector< std::vector< int > > &indices)
Definition JetConstituentsUtils.cc:45
rv::RVec< FCCAnalysesJetConstituentsData > get_dndx(const rv::RVec< FCCAnalysesJetConstituents > &jetConstituents, const TrackUtils::TrackDqdxHandler &dNdxHandler, const rv::RVec< edm4hep::TrackData > &trackColl, const rv::RVec< FCCAnalysesJetConstituentsData > isJetConstChargedHad)
Obtain dNdx corresponding to the provided jet constituents.
Definition JetConstituentsUtils.cc:360
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > jets_sorting_on_nconst(const rv::RVec< edm4hep::ReconstructedParticleData > &)
Definition JetConstituentsUtils.cc:162
rv::RVec< FCCAnalysesJetConstituentsData > get_Sip3dSig(const rv::RVec< FCCAnalysesJetConstituentsData > &Sip3dVals, const rv::RVec< FCCAnalysesJetConstituentsData > &err2_D0, const rv::RVec< FCCAnalysesJetConstituentsData > &err2_Z0)
Definition JetConstituentsUtils.cc:597
rv::RVec< FCCAnalysesJetConstituentsData > XPtoPar_dz(const rv::RVec< FCCAnalysesJetConstituents > &, const ROOT::VecOps::RVec< edm4hep::TrackState > &, const TLorentzVector &V, const float &)
Definition JetConstituentsUtils.cc:234
rv::RVec< FCCAnalysesJetConstituentsData > get_PIDs_cluster(const ROOT::VecOps::RVec< int > recin, const ROOT::VecOps::RVec< int > mcin, const rv::RVec< edm4hep::ReconstructedParticleData > &RecPart, const rv::RVec< edm4hep::MCParticleData > &Particle, const std::vector< std::vector< int > > &indices)
Definition JetConstituentsUtils.cc:1066
rv::RVec< FCCAnalysesJetConstituentsData > get_phirel_cluster(const rv::RVec< fastjet::PseudoJet > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs)
Definition JetConstituentsUtils.cc:1017
rv::RVec< FCCAnalysesJetConstituentsData > get_Sip2dVal(const rv::RVec< edm4hep::ReconstructedParticleData > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:397
rv::RVec< FCCAnalysesJetConstituentsData > get_erel_log_cluster(const rv::RVec< fastjet::PseudoJet > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs)
Definition JetConstituentsUtils.cc:882
rv::RVec< FCCAnalysesJetConstituentsData > get_omega_phi0_cov(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:342
rv::RVec< FCCAnalysesJetConstituentsData > get_d0_z0_cov(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:300
rv::RVec< FCCAnalysesJetConstituentsData > get_isEl(const rv::RVec< FCCAnalysesJetConstituents > &jcs)
Definition JetConstituentsUtils.cc:1088
rv::RVec< double > compute_residue_py(const rv::RVec< TLorentzVector > &tlv_jet, const rv::RVec< TLorentzVector > &sum_tlv_jcs)
Definition JetConstituentsUtils.cc:1338
rv::RVec< FCCAnalysesJetConstituentsData > get_tanlambda_phi0_cov(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:318
rv::RVec< FCCAnalysesJetConstituentsData > get_Sip2dSig(const rv::RVec< FCCAnalysesJetConstituentsData > &Sip2dVals, const rv::RVec< FCCAnalysesJetConstituentsData > &err2_D0)
The functions get_Sip2dSig and get_Sip2dVal can be made independent; I passed to the former the resul...
Definition JetConstituentsUtils.cc:485
rv::RVec< FCCAnalysesJetConstituentsData > get_erel_log(const rv::RVec< edm4hep::ReconstructedParticleData > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs)
Definition JetConstituentsUtils.cc:863
rv::RVec< FCCAnalysesJetConstituentsData > get_Sip3dVal(const rv::RVec< edm4hep::ReconstructedParticleData > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:508
rv::RVec< FCCAnalysesJetConstituentsData > get_phi0(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:207
rv::RVec< FCCAnalysesJetConstituentsData > get_tanlambda_d0_cov(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:324
rv::RVec< int > count_type(const rv::RVec< FCCAnalysesJetConstituentsData > &isType)
Definition JetConstituentsUtils.cc:1226
rv::RVec< double > compute_residue_phi(const rv::RVec< TLorentzVector > &tlv_jet, const rv::RVec< TLorentzVector > &sum_tlv_jcs)
Definition JetConstituentsUtils.cc:1373
rv::RVec< FCCAnalysesJetConstituentsData > get_JetDistVal_cluster(const rv::RVec< fastjet::PseudoJet > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:654
rv::RVec< FCCAnalysesJetConstituentsData > get_omega_d0_cov(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:348
rv::RVec< FCCAnalysesJetConstituentsData > get_phi0_z0_cov(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:312
rv::RVec< FCCAnalysesJetConstituentsData > get_JetDistSig(const rv::RVec< FCCAnalysesJetConstituentsData > &JetDistVal, const rv::RVec< FCCAnalysesJetConstituentsData > &err2_D0, const rv::RVec< FCCAnalysesJetConstituentsData > &err2_Z0)
Definition JetConstituentsUtils.cc:721
rv::RVec< FCCAnalysesJetConstituentsData > XPtoPar_dxy(const rv::RVec< FCCAnalysesJetConstituents > &, const ROOT::VecOps::RVec< edm4hep::TrackState > &, const TLorentzVector &V, const float &)
Definition JetConstituentsUtils.cc:225
rv::RVec< int > count_consts(rv::RVec< FCCAnalysesJetConstituents > jets)
Definition JetConstituentsUtils.cc:1216
rv::RVec< float > FCCAnalysesJetConstituentsData
Definition JetConstituentsUtils.h:35
rv::RVec< FCCAnalysesJetConstituentsData > get_tanlambda_cov(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:294
rv::RVec< FCCAnalysesJetConstituentsData > get_d0_cov(const rv::RVec< FCCAnalysesJetConstituents > &, const ROOT::VecOps::RVec< edm4hep::TrackState > &)
Definition JetConstituentsUtils.cc:276
rv::RVec< FCCAnalysesJetConstituentsData > get_Sip2dVal_cluster(const rv::RVec< fastjet::PseudoJet > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:426
rv::RVec< FCCAnalysesJetConstituentsData > get_thetarel_cluster(const rv::RVec< fastjet::PseudoJet > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs)
Definition JetConstituentsUtils.cc:965
rv::RVec< FCCAnalysesJetConstituentsData > get_d0(const rv::RVec< FCCAnalysesJetConstituents > &, const ROOT::VecOps::RVec< edm4hep::TrackState > &)
Definition JetConstituentsUtils.cc:195
rv::RVec< FCCAnalysesJetConstituentsData > get_omega_z0_cov(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:354
rv::RVec< FCCAnalysesJetConstituentsData > get_tanLambda(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:219
rv::RVec< FCCAnalysesJetConstituentsData > get_omega_cov(const rv::RVec< FCCAnalysesJetConstituents > &, const ROOT::VecOps::RVec< edm4hep::TrackState > &)
Definition JetConstituentsUtils.cc:270
rv::RVec< FCCAnalysesJetConstituentsData > get_Sip3dVal_clusterV(const rv::RVec< fastjet::PseudoJet > &jets, const rv::RVec< FCCAnalysesJetConstituentsData > &D0, const rv::RVec< FCCAnalysesJetConstituentsData > &Z0, const rv::RVec< FCCAnalysesJetConstituentsData > &phi0, const float Bz)
Definition JetConstituentsUtils.cc:568
rv::RVec< double > compute_residue_pz(const rv::RVec< TLorentzVector > &tlv_jet, const rv::RVec< TLorentzVector > &sum_tlv_jcs)
Definition JetConstituentsUtils.cc:1349
rv::RVec< FCCAnalysesJetConstituentsData > XPtoPar_C(const rv::RVec< FCCAnalysesJetConstituents > &, const ROOT::VecOps::RVec< edm4hep::TrackState > &, const float &)
Definition JetConstituentsUtils.cc:252
rv::RVec< FCCAnalysesJetConstituentsData > get_Bz(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:120
rv::RVec< FCCAnalysesJetConstituentsData > get_e(const rv::RVec< FCCAnalysesJetConstituents > &)
Definition JetConstituentsUtils.cc:136
rv::RVec< FCCAnalysesJetConstituentsData > get_isChargedHad(const rv::RVec< FCCAnalysesJetConstituents > &jcs)
Definition JetConstituentsUtils.cc:1136
rv::RVec< double > compute_residue_energy(const rv::RVec< TLorentzVector > &tlv_jet, const rv::RVec< TLorentzVector > &sum_tlv_jcs)
Definition JetConstituentsUtils.cc:1316
rv::RVec< double > compute_residue_theta(const rv::RVec< TLorentzVector > &tlv_jet, const rv::RVec< TLorentzVector > &sum_tlv_jcs)
Definition JetConstituentsUtils.cc:1386
rv::RVec< FCCAnalysesJetConstituents > get_constituents(const rv::RVec< FCCAnalysesJetConstituents > &, const rv::RVec< int > &)
Retrieve the constituents of a collection of indexed jets in event.
Definition JetConstituentsUtils.cc:68
rv::RVec< FCCAnalysesJetConstituentsData > get_erel(const rv::RVec< edm4hep::ReconstructedParticleData > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs)
Definition JetConstituentsUtils.cc:901
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > jets_sorting_on_energy(const rv::RVec< edm4hep::ReconstructedParticleData > &)
Definition JetConstituentsUtils.cc:178
int count_jets(rv::RVec< FCCAnalysesJetConstituents > jets)
Definition JetConstituentsUtils.cc:1211
rv::RVec< double > compute_residue_px(const rv::RVec< TLorentzVector > &tlv_jet, const rv::RVec< TLorentzVector > &sum_tlv_jcs)
Definition JetConstituentsUtils.cc:1327
rv::RVec< FCCAnalysesJetConstituentsData > get_Sip3dVal_cluster(const rv::RVec< fastjet::PseudoJet > &jets, const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:538
rv::RVec< FCCAnalysesJetConstituentsData > get_z0_cov(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:282
rv::RVec< FCCAnalysesJetConstituentsData > get_omega_tanlambda_cov(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const ROOT::VecOps::RVec< edm4hep::TrackState > &tracks)
Definition JetConstituentsUtils.cc:336
rv::RVec< FCCAnalysesJetConstituentsData > XPtoPar_phi(const rv::RVec< FCCAnalysesJetConstituents > &, const ROOT::VecOps::RVec< edm4hep::TrackState > &, const TLorentzVector &V, const float &)
Definition JetConstituentsUtils.cc:243
rv::RVec< FCCAnalysesJetConstituentsData > get_charge(const rv::RVec< FCCAnalysesJetConstituents > &)
Definition JetConstituentsUtils.cc:156
FCCAnalysesJetConstituents get_jet_constituents(const rv::RVec< FCCAnalysesJetConstituents > &, int)
Retrieve the constituents of an indexed jet in event.
Definition JetConstituentsUtils.cc:61
rv::RVec< FCCAnalysesJetConstituentsData > get_p(const rv::RVec< FCCAnalysesJetConstituents > &)
Definition JetConstituentsUtils.cc:131
rv::RVec< FCCAnalysesJetConstituentsData > get_type(const rv::RVec< FCCAnalysesJetConstituents > &)
Definition JetConstituentsUtils.cc:151
rv::RVec< FCCAnalysesJetConstituentsData > get_phi(const rv::RVec< FCCAnalysesJetConstituents > &)
Definition JetConstituentsUtils.cc:146
rv::RVec< FCCAnalysesJetConstituentsData > get_mtof(const rv::RVec< FCCAnalysesJetConstituents > &jcs, const rv::RVec< float > &track_L, const rv::RVec< edm4hep::TrackData > &trackdata, const rv::RVec< edm4hep::TrackerHit3DData > &trackerhits, const rv::RVec< edm4hep::ClusterData > &gammadata, const rv::RVec< edm4hep::ClusterData > &nhdata, const rv::RVec< edm4hep::CalorimeterHitData > &calohits, const TLorentzVector &V)
Definition JetConstituentsUtils.cc:759
rv::RVec< FCCAnalysesJetConstituentsData > XPtoPar_ct(const rv::RVec< FCCAnalysesJetConstituents > &, const ROOT::VecOps::RVec< edm4hep::TrackState > &, const float &)
Definition JetConstituentsUtils.cc:260
rv::RVec< FCCAnalysesJetConstituentsData > get_PIDs(const ROOT::VecOps::RVec< int > recin, const ROOT::VecOps::RVec< int > mcin, const rv::RVec< edm4hep::ReconstructedParticleData > &RecPart, const rv::RVec< edm4hep::MCParticleData > &Particle, const rv::RVec< edm4hep::ReconstructedParticleData > &Jets)
Definition JetConstituentsUtils.cc:1045
rv::RVec< double > compute_residue_pt(const rv::RVec< TLorentzVector > &tlv_jet, const rv::RVec< TLorentzVector > &sum_tlv_jcs)
Definition JetConstituentsUtils.cc:1360
rv::RVec< FCCAnalysesJetConstituentsData > get_Sip2dVal_clusterV(const rv::RVec< fastjet::PseudoJet > &jets, const rv::RVec< FCCAnalysesJetConstituentsData > &D0, const rv::RVec< FCCAnalysesJetConstituentsData > &phi0, const float Bz)
Definition JetConstituentsUtils.cc:455
FCC analyzers collection.
Definition Algorithms.h:15
Definition JetConstituentsUtils.h:20
edm4hep::TrackerHitData TrackerHit3DData
Definition JetConstituentsUtils.h:21