ATLAS Offline Software
Loading...
Searching...
No Matches
SubjetMakerTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
7#include "fastjet/ClusterSequence.hh"
8
9SubjetMakerTool::SubjetMakerTool(const std::string& name) :
11{
12 declareProperty("type", m_type = "AntiKt");
13 declareProperty("R", m_R = 0.2);
14 declareProperty("PtCut", m_minPt = 20e3);
15}
16
18 if(!checkForConstituents(jet)) return 1;
19
20 std::vector<float> s_e;
21 std::vector<float> s_px;
22 std::vector<float> s_py;
23 std::vector<float> s_pz;
24
25 std::vector<fastjet::PseudoJet> constit_pseudojets = jet::JetConstituentFiller::constituentPseudoJets(jet);
26 fastjet::PseudoJet pjet = fastjet::join(constit_pseudojets);
27
28 fastjet::JetAlgorithm jetalg = fastjet::antikt_algorithm;
29 if (m_type == "AntiKt") {
30 jetalg = fastjet::antikt_algorithm;
31 } else if (m_type == "Kt") {
32 jetalg = fastjet::kt_algorithm;
33 } else if (m_type == "CamKt") {
34 jetalg = fastjet::cambridge_algorithm;
35 }
36
37 fastjet::JetDefinition microjet_def(jetalg, m_R);
38 fastjet::ClusterSequence microjet_cs(pjet.constituents(), microjet_def);
39 std::vector<fastjet::PseudoJet> microjets = fastjet::sorted_by_pt(microjet_cs.inclusive_jets(m_minPt));
40 for (size_t z = 0; z < microjets.size(); ++z) {
41 s_e.push_back(microjets[z].e());
42 s_px.push_back(microjets[z].px());
43 s_py.push_back(microjets[z].py());
44 s_pz.push_back(microjets[z].pz());
45 }
46
47 std::stringstream ss;
48 ss.str("");
49 ss << "Subjet" << m_type << (int) (m_R*100) << "_e";
50 jet.setAttribute(ss.str().c_str(), s_e);
51 ss.str("");
52 ss << "Subjet" << m_type << (int) (m_R*100) << "_px";
53 jet.setAttribute(ss.str().c_str(), s_px);
54 ss.str("");
55 ss << "Subjet" << m_type << (int) (m_R*100) << "_py";
56 jet.setAttribute(ss.str().c_str(), s_py);
57 ss.str("");
58 ss << "Subjet" << m_type << (int) (m_R*100) << "_pz";
59 jet.setAttribute(ss.str().c_str(), s_pz);
60
61 return 0;
62}
static Double_t ss
#define z
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
JetSubStructureMomentToolsBase(const std::string &name)
bool checkForConstituents(const xAOD::Jet &jet) const
std::string m_type
int modifyJet(xAOD::Jet &jet) const
Modify a single jet. This is obsolete and set to be removed.
SubjetMakerTool(const std::string &name)
static PseudoJetVector constituentPseudoJets(const xAOD::Jet &jet, bool ignoreGhosts=true, bool requireJetStructure=false)
Returns the jet's constituents as a vector of PseudoJet if ignoreGhosts==true, ghost constituents are...
Jet_v1 Jet
Definition of the current "jet version".