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

#include <Qw.h>

Inheritance diagram for JetSubStructureUtils::Qw:
Collaboration diagram for JetSubStructureUtils::Qw:

Public Types

enum  QwVariant { Normal , MassCut , SmallSubjets }
 describes the differnt calculation variants More...

Public Member Functions

 Qw (QwVariant mode=Normal, double p=-1, fastjet::JetAlgorithm jetalg=fastjet::kt_algorithm)
virtual double result (const fastjet::PseudoJet &jet) const
 SubstructureCalculator interface.
virtual double result (std::vector< fastjet::PseudoJet > &constituents) const
 Implementation of the calculation.
virtual TOut result (const xAOD::Jet &jet) const

Private Attributes

QwVariant m_vMode
fastjet::JetDefinition m_jetdef
double m_massCut

Detailed Description

Definition at line 34 of file Qw.h.

Member Enumeration Documentation

◆ QwVariant

describes the differnt calculation variants

Enumerator
Normal 
MassCut 
SmallSubjets 

Definition at line 38 of file Qw.h.

Constructor & Destructor Documentation

◆ Qw()

Qw::Qw ( QwVariant mode = Normal,
double p = -1,
fastjet::JetAlgorithm jetalg = fastjet::kt_algorithm )

Definition at line 13 of file Qw.cxx.

13 : m_vMode(mode), m_jetdef(), m_massCut(-1) {
14
15 switch (mode ) {
16 case Normal:
17 case MassCut:
18 m_jetdef = fastjet::JetDefinition(jetalg, 1.5) ;
19 m_massCut = p;
20 break;
21 case SmallSubjets:
22 if(p==-1) p = 0.2;
23 m_jetdef = fastjet::JetDefinition(jetalg, p) ;
24 break;
25 }
26}
QwVariant m_vMode
Definition Qw.h:55
fastjet::JetDefinition m_jetdef
Definition Qw.h:56

Member Function Documentation

◆ result() [1/3]

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

SubstructureCalculator interface.

Definition at line 28 of file Qw.cxx.

29{
30 std::vector<fastjet::PseudoJet> constituents = jet.constituents();
31 return result(constituents);
32}
virtual double result(const fastjet::PseudoJet &jet) const
SubstructureCalculator interface.
Definition Qw.cxx:28

◆ result() [2/3]

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 }
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

◆ result() [3/3]

double Qw::result ( std::vector< fastjet::PseudoJet > & constituents) const
virtual

Implementation of the calculation.

Definition at line 53 of file Qw.cxx.

53 {
54
55 // deal with extreme cases first.
56 // early exit in this cases.
57 size_t nconst = constituents.size();
58 if( nconst < 3 ) return 0;
59
60 // Build the subjets
61 fastjet::ClusterSequence cs(constituents, m_jetdef);
62
63 double scaleF = 1.;
64
65 std::vector<fastjet::PseudoJet> outjets;
66
67 switch( m_vMode) {
68 case Normal :
69 outjets= cs.exclusive_jets(3);
70 //std::cout << " Normal mode : scaleF " << scaleF << " m1="<< outjets[0].m() << std::endl;
71
72 break;
73
74 case MassCut: {
75 outjets= cs.exclusive_jets(3);
76 fastjet::PseudoJet sumTot(0,0,0,0) ;
77 sumTot = outjets[0] + outjets[1] + outjets[2];
78
79 // make sure we deal only with subjets with mass < m_massCut
80 std::vector<fastjet::PseudoJet> tmp_subjets; tmp_subjets.reserve(4);
81 decomposeToJetMass(outjets[0], m_massCut , tmp_subjets);
82 decomposeToJetMass(outjets[1], m_massCut , tmp_subjets);
83 decomposeToJetMass(outjets[2], m_massCut , tmp_subjets);
84 // Take the 3 highest pt subjets.
85 outjets = fastjet::sorted_by_pt( tmp_subjets) ;
86 outjets.resize(3);
87
88
89 // resum the 3 highest pt subjets to compute the scale factor.
90 fastjet::PseudoJet sumSubJ(0,0,0,0);
91 sumSubJ = outjets[0] + outjets[1] + outjets[2];
92 scaleF = sumTot.pt()/sumSubJ.pt();
93 // std::cout << " MassCut mode : scaleF " << scaleF << " m1="<< outjets[0].m() << std::endl;
94
95 }
96 break;
97 case SmallSubjets: {
98 outjets = fastjet::sorted_by_pt(cs.inclusive_jets());
99
100 // deal with extreme case first.
101 // early exit in this case.
102 size_t nsubjets = outjets.size();
103 if( nsubjets < 3 ) return 0;
104
105 // compute the scale factor : as we choose only the 1st 3 subjets,
106 // we scale back the contribution of these 3 subjets to the total Pt.
107 fastjet::PseudoJet sumTot(0,0,0,0), sumSubJ(0,0,0,0);
108 for( fastjet::PseudoJet& c: constituents) sumTot+=c;
109
110 // take the 3 highest subjets only
111 outjets.resize(3);
112
113 for( fastjet::PseudoJet& c: outjets ) sumSubJ+=c;
114 scaleF = sumTot.pt()/sumSubJ.pt();
115 // std::cout << " SmallSubjets mode : scaleF " << scaleF << " m1="<< outjets[0].m() << std::endl;
116 }
117 break;
118 }
119
120 double m12 = (outjets[0]+outjets[1]).m();
121 double m23 = (outjets[2]+outjets[1]).m();
122 double m13 = (outjets[2]+outjets[0]).m();
123
124 // std::cout << " ---------> qw="<< scaleF*std::min(m12, std::min(m23,m13) ) << std::endl;
125 return scaleF*std::min(m12, std::min(m23,m13) );
126}

Member Data Documentation

◆ m_jetdef

fastjet::JetDefinition JetSubStructureUtils::Qw::m_jetdef
private

Definition at line 56 of file Qw.h.

◆ m_massCut

double JetSubStructureUtils::Qw::m_massCut
private

Definition at line 58 of file Qw.h.

◆ m_vMode

QwVariant JetSubStructureUtils::Qw::m_vMode
private

Definition at line 55 of file Qw.h.


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