ATLAS Offline Software
Loading...
Searching...
No Matches
JetSubStructureUtils::SubjetFinder Class Reference

#include <SubjetFinder.h>

Inheritance diagram for JetSubStructureUtils::SubjetFinder:
Collaboration diagram for JetSubStructureUtils::SubjetFinder:

Public Member Functions

 SubjetFinder (fastjet::JetAlgorithm fj_jetalg=fastjet::kt_algorithm, float jet_radius=0.3, float pt_min=5000, int exclusive_njets=-1)
virtual std::vector< fastjet::PseudoJet > result (const fastjet::PseudoJet &jet) const
virtual TOut result (const xAOD::Jet &jet) const

Private Attributes

fastjet::JetAlgorithm m_fj_jetalg
float m_jetrad
float m_ptmin
int m_exclusivenjets

Detailed Description

Definition at line 13 of file SubjetFinder.h.

Constructor & Destructor Documentation

◆ SubjetFinder()

SubjetFinder::SubjetFinder ( fastjet::JetAlgorithm fj_jetalg = fastjet::kt_algorithm,
float jet_radius = 0.3,
float pt_min = 5000,
int exclusive_njets = -1 )

Definition at line 15 of file SubjetFinder.cxx.

15 :
16 m_fj_jetalg(fj_jetalg), m_jetrad(jet_radius), m_ptmin(pt_min), m_exclusivenjets(exclusive_njets)
17{
18}
fastjet::JetAlgorithm m_fj_jetalg

Member Function Documentation

◆ result() [1/2]

vector< fastjet::PseudoJet > SubjetFinder::result ( const fastjet::PseudoJet & jet) const
virtual

Definition at line 20 of file SubjetFinder.cxx.

21{
22 vector<fastjet::PseudoJet> constit_pseudojets = jet.constituents();
23 vector<fastjet::PseudoJet> subjets;
24 if(constit_pseudojets.empty()) {
25 cout << "Warning in SubjetFinder: jet has no constituents" << endl;
26 return subjets;
27 }
28
29 fastjet::ClusterSequence *clust_seq = nullptr;
30
31 if (m_fj_jetalg==fastjet::ee_kt_algorithm) {
32 fastjet::JetDefinition jet_def = fastjet::JetDefinition(fastjet::ee_kt_algorithm);
33 clust_seq = new fastjet::ClusterSequence(constit_pseudojets, jet_def);
34 }else if (m_fj_jetalg==fastjet::plugin_algorithm) {
35 fastjet::JetDefinition jet_def = fastjet::JetDefinition(new fastjet::EECambridgePlugin(m_jetrad));
36 clust_seq = new fastjet::ClusterSequence(constit_pseudojets, jet_def);
37 }else{
38 fastjet::JetDefinition jet_def = fastjet::JetDefinition(m_fj_jetalg, m_jetrad, fastjet::E_scheme, fastjet::Best);
39 clust_seq = new fastjet::ClusterSequence(constit_pseudojets, jet_def);
40 }
41
42 if(m_exclusivenjets < 0) { // Inclusive
43 subjets = fastjet::sorted_by_pt(clust_seq->inclusive_jets(m_ptmin));
44 }
45 else {
46 subjets = fastjet::sorted_by_pt(clust_seq->exclusive_jets_up_to(m_exclusivenjets));
47 }
48
49 if(subjets.empty()) {
50 delete clust_seq;
51 }
52 else {
53 clust_seq->delete_self_when_unused();
54 }
55
56 return subjets;
57}

◆ result() [2/2]

virtual TOut JetSubStructureUtils::SubstructureCalculator< TOut >::result ( const xAOD::Jet & jet) const
inlinevirtual

Reimplemented from JetSubStructureUtils::SubstructureCalculator< std::vector< fastjet::PseudoJet > >.

Definition at line 25 of file SubstructureCalculator.h.

25 {
26 // PS 4/18 master developent
27 // std::vector<fastjet::PseudoJet> constit_pseudojets =
28 // jet::JetConstituentFiller::constituentPseudoJets(jet);
29
30 std::vector<fastjet::PseudoJet> constit_pseudojets;
31 std::transform(jet.getConstituents().begin(),
32 jet.getConstituents().end(),
33 std::back_inserter(constit_pseudojets),
34 [](const auto& c){
35 const xAOD::IParticle* ip = c->rawConstituent();
36 return
37 // fastjet::PseudoJet((c->rawConstituent())->p4());
38 fastjet::PseudoJet(ip->p4());
39 });
40
41 fastjet::PseudoJet pjet = fastjet::join(constit_pseudojets);
42
43 return result(pjet);
44 }
virtual std::vector< fastjet::PseudoJet > result(const fastjet::PseudoJet &jet) const
iterator begin() const
iterator on the first constituent
iterator end() const
iterator after the last constituent
JetConstituentVector getConstituents() const
Return a vector of consituents. The object behaves like vector<const IParticle*>. See JetConstituentV...
Definition Jet_v1.cxx:147

Member Data Documentation

◆ m_exclusivenjets

int JetSubStructureUtils::SubjetFinder::m_exclusivenjets
private

Definition at line 24 of file SubjetFinder.h.

◆ m_fj_jetalg

fastjet::JetAlgorithm JetSubStructureUtils::SubjetFinder::m_fj_jetalg
private

Definition at line 21 of file SubjetFinder.h.

◆ m_jetrad

float JetSubStructureUtils::SubjetFinder::m_jetrad
private

Definition at line 22 of file SubjetFinder.h.

◆ m_ptmin

float JetSubStructureUtils::SubjetFinder::m_ptmin
private

Definition at line 23 of file SubjetFinder.h.


The documentation for this class was generated from the following files: