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

#include <Dipolarity.h>

Inheritance diagram for JetSubStructureUtils::Dipolarity:
Collaboration diagram for JetSubStructureUtils::Dipolarity:

Public Member Functions

 Dipolarity (unsigned int n1, unsigned int n2, bool exclusive, double subjetR)
virtual double result (const fastjet::PseudoJet &jet) const
virtual TOut result (const xAOD::Jet &jet) const

Private Member Functions

double dipolarity (std::vector< fastjet::PseudoJet > &constit_pseudojets, const fastjet::PseudoJet &jet1, const fastjet::PseudoJet &jet2) const

Private Attributes

bool m_exclusive
unsigned int m_n1
unsigned int m_n2
double m_subjetR

Detailed Description

Definition at line 12 of file Dipolarity.h.

Constructor & Destructor Documentation

◆ Dipolarity()

Dipolarity::Dipolarity ( unsigned int n1,
unsigned int n2,
bool exclusive,
double subjetR )

Definition at line 15 of file Dipolarity.cxx.

15 :
16 m_exclusive(exclusive),
17 m_n1(n1),
18 m_n2(n2),
19 m_subjetR(subjetR) {
20
21 if(m_n2 <= m_n1 || m_n1 == 0) {
22 cout << "dipolarity was called with invalid arguments" << endl;
23 }
24}

Member Function Documentation

◆ dipolarity()

double Dipolarity::dipolarity ( std::vector< fastjet::PseudoJet > & constit_pseudojets,
const fastjet::PseudoJet & jet1,
const fastjet::PseudoJet & jet2 ) const
private

Definition at line 62 of file Dipolarity.cxx.

64{
65 float dipolarity = 0;
66 float sumpt = 0;
67
68 TVector2 v12(jet2.eta() - jet1.eta(), jet1.delta_phi_to(jet2));
69 if(v12.Mod2() < 0.001) return -1;
70
71 for(unsigned int iConstit = 0; iConstit < constit_pseudojets.size(); iConstit++) {
72 fastjet::PseudoJet constituent = constit_pseudojets[iConstit];
73
74 sumpt += constituent.perp();
75
76 TVector2 v;
77
78 v.Set(constituent.eta() - jet1.eta(), jet1.delta_phi_to(constituent));
79 float test = v * v12.Unit();
80
81 if(test < 0)
82 dipolarity += constituent.perp()*v.Mod2();
83 else if(test < v12.Mod())
84 dipolarity += constituent.perp()*(v.Mod2() - pow(v * v12.Unit(), 2));
85 else {
86 v.Set(constituent.eta() - jet2.eta(), jet2.delta_phi_to(constituent));
87 dipolarity += constituent.perp()*v.Mod2();
88 }
89 }
90
91 if(sumpt < 0.001) return -1;
92 return dipolarity/(sumpt*v12.Mod2());
93}
constexpr int pow(int base, int exp) noexcept
double dipolarity(std::vector< fastjet::PseudoJet > &constit_pseudojets, const fastjet::PseudoJet &jet1, const fastjet::PseudoJet &jet2) const

◆ result() [1/2]

double Dipolarity::result ( const fastjet::PseudoJet & jet) const
virtual

Definition at line 26 of file Dipolarity.cxx.

27{
28 vector<fastjet::PseudoJet> constit_pseudojets = jet.constituents();
29 if(constit_pseudojets.empty()) return -1;
30
31 vector<fastjet::PseudoJet> kt_subjets;
32 fastjet::PseudoJet jet1, jet2;
33
34 if(!m_exclusive) {
35 fastjet::JetDefinition jet_def = fastjet::JetDefinition(fastjet::kt_algorithm, m_subjetR,
36 fastjet::E_scheme, fastjet::Best);
37 fastjet::ClusterSequence kt_clust_seq(constit_pseudojets, jet_def);
38 kt_subjets = fastjet::sorted_by_pt(kt_clust_seq.inclusive_jets(5000.0));
39
40 if(kt_subjets.size() < m_n2) {
41 //Not enough subjets
42 return 0.0;
43 }
44 }
45 else {
46 if(m_n2 > constit_pseudojets.size()) {
47 //We were asked to calculate exclusiveDipolarity, but there are not enough constituents
48 return 0.0;
49 }
50
51 fastjet::JetDefinition jet_def = fastjet::JetDefinition(fastjet::kt_algorithm, 1.5,
52 fastjet::E_scheme, fastjet::Best);
53 fastjet::ClusterSequence kt_clust_seq(constit_pseudojets, jet_def);
54 kt_subjets = fastjet::sorted_by_pt(kt_clust_seq.exclusive_jets((int)m_n2));
55 }
56
57 jet1 = kt_subjets[m_n1-1];
58 jet2 = kt_subjets[m_n2-1];
59 return dipolarity(constit_pseudojets, jet1, jet2);
60}

◆ result() [2/2]

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

Reimplemented from JetSubStructureUtils::SubstructureCalculator< double >.

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 double 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_exclusive

bool JetSubStructureUtils::Dipolarity::m_exclusive
private

Definition at line 23 of file Dipolarity.h.

◆ m_n1

unsigned int JetSubStructureUtils::Dipolarity::m_n1
private

Definition at line 24 of file Dipolarity.h.

◆ m_n2

unsigned int JetSubStructureUtils::Dipolarity::m_n2
private

Definition at line 25 of file Dipolarity.h.

◆ m_subjetR

double JetSubStructureUtils::Dipolarity::m_subjetR
private

Definition at line 26 of file Dipolarity.h.


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