Skip to content

Commit

Permalink
Merge pull request #107 from gavinsalam/master
Browse files Browse the repository at this point in the history
Adding JetClusteringUtils::set_pseudoJets_xyzm and illustrating usage in examples/FCCee/top/hadronic/analysis.py
  • Loading branch information
clementhelsens authored Sep 20, 2021
2 parents 1703add + edcabba commit 3122187
Show file tree
Hide file tree
Showing 3 changed files with 36 additions and 3 deletions.
19 changes: 19 additions & 0 deletions analyzers/dataframe/JetClusteringUtils.cc
Original file line number Diff line number Diff line change
Expand Up @@ -25,6 +25,25 @@ std::vector<fastjet::PseudoJet> JetClusteringUtils::set_pseudoJets(ROOT::VecOps:
return result;
}

std::vector<fastjet::PseudoJet> JetClusteringUtils::set_pseudoJets_xyzm(ROOT::VecOps::RVec<float> px,
ROOT::VecOps::RVec<float> py,
ROOT::VecOps::RVec<float> pz,
ROOT::VecOps::RVec<float> m) {
std::vector<fastjet::PseudoJet> result;
unsigned index = 0;
for (size_t i = 0; i < px.size(); ++i) {
double px_d = px.at(i);
double py_d = py.at(i);
double pz_d = pz.at(i);
double m_d = m.at(i);
double E_d = sqrt(px_d*px_d + py_d*py_d + pz_d*pz_d + m_d*m_d);
result.emplace_back(px_d, py_d, pz_d, E_d);
result.back().set_user_index(index);
++index;
}
return result;
}


ROOT::VecOps::RVec<float> JetClusteringUtils::get_px(ROOT::VecOps::RVec<fastjet::PseudoJet> in){
ROOT::VecOps::RVec<float> result;
Expand Down
13 changes: 13 additions & 0 deletions analyzers/dataframe/JetClusteringUtils.h
Original file line number Diff line number Diff line change
Expand Up @@ -33,6 +33,19 @@ namespace JetClusteringUtils{
ROOT::VecOps::RVec<float> pz,
ROOT::VecOps::RVec<float> e);

/** Set fastjet pseudoJet for later reconstruction using px, py, pz and m
*
* This version is to be preferred over the px,py,pz,E version when m is known
* accurately, because it uses double precision to reconstruct the energy,
* reducing the size of rounding errors on FastJet calculations (e.g. of
* PseudoJet masses)
*
*/
std::vector<fastjet::PseudoJet> set_pseudoJets_xyzm(ROOT::VecOps::RVec<float> px,
ROOT::VecOps::RVec<float> py,
ROOT::VecOps::RVec<float> pz,
ROOT::VecOps::RVec<float> m);

/** Get fastjet pseudoJet after reconstruction from FCCAnalyses jets*/
ROOT::VecOps::RVec<fastjet::PseudoJet> get_pseudoJets(FCCAnalysesJet);

Expand Down
7 changes: 4 additions & 3 deletions examples/FCCee/top/hadronic/analysis.py
Original file line number Diff line number Diff line change
Expand Up @@ -36,10 +36,11 @@ def run(self):
.Define("RP_px", "ReconstructedParticle::get_px(ReconstructedParticles)")
.Define("RP_py", "ReconstructedParticle::get_py(ReconstructedParticles)")
.Define("RP_pz", "ReconstructedParticle::get_pz(ReconstructedParticles)")
.Define("RP_e", "ReconstructedParticle::get_e(ReconstructedParticles)")
.Define("RP_m", "ReconstructedParticle::get_m(ReconstructedParticles)")

#build pseudo jets with the RP
.Define("pseudo_jets", "JetClusteringUtils::set_pseudoJets(RP_px, RP_py, RP_pz, RP_e)")
#build pseudo jets with the RP, using the interface that takes px,py,pz,m for better
#handling of rounding errors
.Define("pseudo_jets", "JetClusteringUtils::set_pseudoJets_xyzm(RP_px, RP_py, RP_pz, RP_m)")

#run jet clustering with all reconstructed particles. kt_algorithm, R=0.5, exclusive clustering, exactly 4 jets, E0-scheme
.Define("FCCAnalysesJets_kt", "JetClustering::clustering_kt(0.5, 2, 4, 0, 10)(pseudo_jets)")
Expand Down

0 comments on commit 3122187

Please sign in to comment.