ATLAS Offline Software
Loading...
Searching...
No Matches
DerivationFramework::RCJetSubstructureAug Class Reference

#include <RCJetSubstructureAug.h>

Inheritance diagram for DerivationFramework::RCJetSubstructureAug:
Collaboration diagram for DerivationFramework::RCJetSubstructureAug:

Public Member Functions

virtual StatusCode initialize () override final
virtual StatusCode addBranches (const EventContext &ctx) const override final

Private Types

using WDHK = SG::WriteDecorHandleKey<xAOD::JetContainer>

Private Attributes

StringProperty m_streamName { this, "StreamName", "", "Name of the stream" }
Gaudi::Property< std::vector< std::string > > m_ghostNames { this, "GhostConstitNames", {"GhostLCTopo"}, "Names of the ghost constituents for substructure computation"}
StringProperty m_selectionString { this, "SelectionString", "", "Selection to apply to jet"}
StringProperty m_grooming { this, "Grooming", "", "Name of gromming technic to apply (Trimming or SofDrop)"}
Gaudi::Property< float > m_rclus {this, "RClusTrim", 0.3 , "R for reclustering (0 for none)"}
Gaudi::Property< float > m_ptfrac {this, "PtFracTrim", 0.03, "pT fraction for retaining subjets"}
Gaudi::Property< float > m_beta {this, "BetaSoft", 1. , "How much to consider angular dependence"}
Gaudi::Property< float > m_zcut {this, "ZcutSoft", 1. , "pT fraction for retaining subjets"}
Gaudi::Property< float > m_R0 {this, "R0Soft", 1. , "Normalization of angular distance, usually the characteristic jet radius (default R0 = 1)"}
SG::ReadHandleKey< xAOD::JetContainerm_jetKey { this, "JetContainerKey", ""}
WDHK m_dec_Qw { this, "dec_Qw", m_jetKey, "Qw_" }
 Qw decorator.
WDHK m_dec_Tau1 { this, "dec_Tau1", m_jetKey, "Tau1_" }
 Nsubjetiness decorators.
WDHK m_dec_Tau2 { this, "dec_Tau2", m_jetKey, "Tau2_" }
WDHK m_dec_Tau3 { this, "dec_Tau3", m_jetKey, "Tau3_" }
WDHK m_dec_Tau4 { this, "dec_Tau4", m_jetKey, "Tau4_" }
WDHK m_dec_Tau21 { this, "dec_Tau21", m_jetKey, "Tau21_" }
WDHK m_dec_Tau32 { this, "dec_Tau32", m_jetKey, "Tau32_" }
WDHK m_dec_Split12 { this, "dec_Split12", m_jetKey, "Split12_" }
 KtSplittingScale decorators.
WDHK m_dec_Split23 { this, "dec_Split23", m_jetKey, "Split23_" }
WDHK m_dec_Split34 { this, "dec_Split34", m_jetKey, "Split34_" }
WDHK m_dec_ECF1 { this, "dec_ECF1", m_jetKey, "ECF1_" }
 Energy correlation factors decorators.
WDHK m_dec_ECF2 { this, "dec_ECF2", m_jetKey, "ECF2_" }
WDHK m_dec_ECF3 { this, "dec_ECF3", m_jetKey, "ECF3_" }
WDHK m_dec_ECF4 { this, "dec_ECF4", m_jetKey, "ECF4_" }
WDHK m_dec_C2 { this, "dec_C2", m_jetKey, "C2_" }
WDHK m_dec_D2 { this, "dec_D2", m_jetKey, "D2_" }
WDHK m_dec_pT { this, "dec_pT", m_jetKey, "pT_" }
 Reclustered jets information decorators.
WDHK m_dec_m { this, "dec_m", m_jetKey, "m_" }
WDHK m_dec_NConstits { this, "dec_NConstits", m_jetKey, "NConstits_" }
WDHK m_dec_eta { this, "dec_eta", m_jetKey, "eta_" }
WDHK m_dec_phi { this, "dec_phi", m_jetKey, "phi_" }
WDHK m_dec_timing { this, "dec_timing", m_jetKey, "timing_" }
 Timing information.
std::optional< fastjet::Filter > m_trimmer
std::optional< fastjet::contrib::SoftDrop > m_softdropper

Detailed Description

Definition at line 28 of file RCJetSubstructureAug.h.

Member Typedef Documentation

◆ WDHK

Member Function Documentation

◆ addBranches()

StatusCode DerivationFramework::RCJetSubstructureAug::addBranches ( const EventContext & ctx) const
finaloverridevirtual

Definition at line 69 of file RCJetSubstructureAug.cxx.

70{
71
72 SG::ReadHandle<xAOD::JetContainer> jets(m_jetKey,ctx);
73 if (!jets.isValid()) {
74 ATH_MSG_ERROR("No jet collection with name " << m_jetKey.key() << " found in StoreGate!");
75 return StatusCode::FAILURE;
76 }
77 unsigned int nJets(jets->size());
78 std::vector<const xAOD::Jet*> jetToCheck; jetToCheck.clear();
79
80 // Execute the text parser if requested
81 if (!m_selectionString.empty()) {
82 std::vector<int> entries = m_parser->evaluateAsVector();
83 unsigned int nEntries = entries.size();
84 // check the sizes are compatible
85 if (nJets != nEntries ) {
86 ATH_MSG_ERROR("Sizes incompatible! Are you sure your selection string used jets??");
87 return StatusCode::FAILURE;
88 } else {
89 // identify which jets to keep for the thinning check
90 for (unsigned int i=0; i<nJets; ++i) if (entries[i]==1) jetToCheck.push_back((*jets)[i]);
91 }
92 } else {
93 for (unsigned int i=0; i<nJets; ++i) jetToCheck.push_back((*jets)[i]);
94 }
95
96 using WDH = SG::WriteDecorHandle<xAOD::JetContainer, float>;
97 WDH dec_Qw (m_dec_Qw, ctx);
98 WDH dec_Tau1 (m_dec_Tau1, ctx);
99 WDH dec_Tau2 (m_dec_Tau2, ctx);
100 WDH dec_Tau3 (m_dec_Tau3, ctx);
101 WDH dec_Tau4 (m_dec_Tau4, ctx);
102 WDH dec_Tau21 (m_dec_Tau21, ctx);
103 WDH dec_Tau32 (m_dec_Tau32, ctx);
104 WDH dec_Split12 (m_dec_Split12, ctx);
105 WDH dec_Split23 (m_dec_Split23, ctx);
106 WDH dec_Split34 (m_dec_Split34, ctx);
107 WDH dec_ECF1 (m_dec_ECF1, ctx);
108 WDH dec_ECF2 (m_dec_ECF2, ctx);
109 WDH dec_ECF3 (m_dec_ECF3, ctx);
110 WDH dec_ECF4 (m_dec_ECF4, ctx);
111 WDH dec_C2 (m_dec_C2, ctx);
112 WDH dec_D2 (m_dec_D2, ctx);
113 WDH dec_pT (m_dec_pT, ctx);
114 WDH dec_m (m_dec_m, ctx);
115 WDH dec_NConstits (m_dec_NConstits, ctx);
116 WDH dec_eta (m_dec_eta, ctx);
117 WDH dec_phi (m_dec_phi, ctx);
118 WDH dec_timing (m_dec_timing, ctx);
119
120 std::vector<const xAOD::IParticle*> mergedGhostConstits;
121 std::vector<const xAOD::IParticle*> ghosts;
122 std::vector<fastjet::PseudoJet> constituents;
123
124 for (const xAOD::Jet* jet : jetToCheck) {
125
126 // Get list of ghost constituents
127 mergedGhostConstits.clear();
128 for (std::string ghostName: m_ghostNames) {
129 ghosts.clear();
130 ATH_CHECK(jet->getAssociatedObjects(ghostName, ghosts));
131 mergedGhostConstits.insert(mergedGhostConstits.end(), ghosts.begin(), ghosts.end());
132 }
133
134 // Construct list of constituent PseudoJets from constituents and compute timing information
135 constituents.clear();
136 double eTot = 0;
137 double time = 0;
138 for (auto constit : mergedGhostConstits){
139 if (constit==nullptr) {
140 ATH_MSG_INFO("Invalid constituent link, skipping ...");
141 continue;
142 }
143 // Filter constituents on negative energy
144 if (constit->e()<0) {ATH_MSG_INFO("################## Negative energy cluster"); continue;}
145 constituents.push_back( fastjet::PseudoJet(constit->p4()) );
146
147 // Timing info of constituent only for calo constituents
148 if (constit->type() == xAOD::Type::CaloCluster) {
149 auto caloConstit = static_cast<const xAOD::CaloCluster*> (constit);
150 double eConstit = caloConstit->e()* caloConstit->e();
151 time += caloConstit->time()* eConstit;
152 eTot += eConstit;
153 }
154 }
155
156 // Fill timing info
157 if (eTot==0) {
158 dec_timing(*jet) = 0;
159 } else {
160 dec_timing(*jet) = time/eTot;
161 }
162
163 // Get number of constituents
164 int nConstits = constituents.size();
165
166 // Run clustering on constituents if not empty
167 fastjet::PseudoJet groomed_jet;
168 if (nConstits!=0) {
169 auto jet_def = fastjet::JetDefinition(fastjet::antikt_algorithm, 1.5);
170 fastjet::ClusterSequence cs(constituents, jet_def);
171 fastjet::PseudoJet recluster_jet = cs.inclusive_jets(0.0).front();
172
173 // Apply grooming to reclustered jet
174 if (m_grooming=="Trimming"){
175 groomed_jet = m_trimmer->result(recluster_jet);
176 } else if (m_grooming=="SoftDrop"){
177 groomed_jet = m_softdropper->result(recluster_jet);
178 } else {
179 ATH_MSG_DEBUG(" No grooming requested or wrong one, will not apply grooming");
180 groomed_jet = recluster_jet;
181 }
182
183 // update nConstit
184 nConstits = groomed_jet.constituents().size();
185 }
186
187 // Fill substructure vars with default values when no constituents in jet
188 if (nConstits==0) {
189 dec_Qw(*jet) = -999;
190
191 dec_Tau1(*jet) = -999;
192 dec_Tau2(*jet) = -999;
193 dec_Tau3(*jet) = -999;
194 dec_Tau4(*jet) = -999;
195
196 dec_Tau21(*jet) = -999;
197 dec_Tau32(*jet) = -999;
198
199 dec_Split12(*jet) = -999;
200 dec_Split23(*jet) = -999;
201 dec_Split34(*jet) = -999;
202
203 dec_ECF1(*jet) = -999;
204 dec_ECF2(*jet) = -999;
205 dec_ECF3(*jet) = -999;
206 dec_ECF4(*jet) = -999;
207
208 dec_C2(*jet) = -999;
209 dec_D2(*jet) = -999;
210
211 dec_pT(*jet) = 0;
212 dec_m(*jet) = 0;
213 dec_NConstits(*jet) = 0;
214 dec_eta(*jet) = -999;
215 dec_phi(*jet) = -999;
216
217 return StatusCode::SUCCESS;
218 }
219
220 // Save reclustered jet info
221 dec_pT(*jet) = groomed_jet.pt();
222 dec_m(*jet) = groomed_jet.m();
223 dec_NConstits(*jet) = nConstits;
224 dec_eta(*jet) = groomed_jet.eta();
225 dec_phi(*jet) = groomed_jet.phi() - M_PI;
226
227 // Qw
228 JetSubStructureUtils::Qw qw;
229 dec_Qw(*jet) = qw.result(groomed_jet);
230
231 // Nsubjetiness
232 fastjet::contrib::WTA_KT_Axes wta_kt_axes;
233 fastjet::contrib::NormalizedCutoffMeasure normalized_measure(1.0, 1.0, 1000000);
234 JetSubStructureUtils::Nsubjettiness Tau1_wta(1, wta_kt_axes, normalized_measure);
235 JetSubStructureUtils::Nsubjettiness Tau2_wta(2, wta_kt_axes, normalized_measure);
236 JetSubStructureUtils::Nsubjettiness Tau3_wta(3, wta_kt_axes, normalized_measure);
237 JetSubStructureUtils::Nsubjettiness Tau4_wta(4, wta_kt_axes, normalized_measure);
238 float tau1_wta = Tau1_wta.result(groomed_jet);
239 float tau2_wta = Tau2_wta.result(groomed_jet);
240 float tau3_wta = Tau3_wta.result(groomed_jet);
241 float tau4_wta = Tau4_wta.result(groomed_jet);
242 float tau21_wta = -999;
243 float tau32_wta = -999;
244
245 if( tau1_wta > 1e-8 ) {
246 tau21_wta = tau2_wta / tau1_wta;
247 }
248 if( tau2_wta > 1e-8 ) {
249 tau32_wta = tau3_wta / tau2_wta;
250 }
251
252 dec_Tau1(*jet) = tau1_wta;
253 dec_Tau2(*jet) = tau2_wta;
254 dec_Tau3(*jet) = tau3_wta;
255 dec_Tau4(*jet) = tau4_wta;
256 dec_Tau21(*jet) = tau21_wta;
257 dec_Tau32(*jet) = tau32_wta;
258
259 // KtSplittingScale
260 JetSubStructureUtils::KtSplittingScale Split12(1);
261 JetSubStructureUtils::KtSplittingScale Split23(2);
262 JetSubStructureUtils::KtSplittingScale Split34(3);
263 float split12 = Split12.result(groomed_jet);
264 float split23 = Split23.result(groomed_jet);
265 float split34 = Split34.result(groomed_jet);
266
267 dec_Split12(*jet) = split12;
268 dec_Split23(*jet) = split23;
269 dec_Split34(*jet) = split34;
270
271 // EnergyCorrelator
272 JetSubStructureUtils::EnergyCorrelator ECF1(1, 1.0, JetSubStructureUtils::EnergyCorrelator::pt_R);
273 JetSubStructureUtils::EnergyCorrelator ECF2(2, 1.0, JetSubStructureUtils::EnergyCorrelator::pt_R);
274 JetSubStructureUtils::EnergyCorrelator ECF3(3, 1.0, JetSubStructureUtils::EnergyCorrelator::pt_R);
275 JetSubStructureUtils::EnergyCorrelator ECF4(4, 1.0, JetSubStructureUtils::EnergyCorrelator::pt_R);
276 float ecf1 = ECF1.result(groomed_jet);
277 float ecf2 = ECF2.result(groomed_jet);
278 float ecf3 = ECF3.result(groomed_jet);
279 float ecf4 = ECF4.result(groomed_jet);
280 float c2 = -999;
281 float d2 = -999;
282
283 if( ecf2 > 1e-8 ) {
284 c2 = ecf3 * ecf1 / pow( ecf2, 2.0 );
285 d2 = ecf3 * pow( ecf1, 3.0 ) / pow( ecf2, 3.0 );
286 }
287
288 dec_ECF1(*jet) = ecf1;
289 dec_ECF2(*jet) = ecf2;
290 dec_ECF3(*jet) = ecf3;
291 dec_ECF4(*jet) = ecf4;
292 dec_C2(*jet) = c2;
293 dec_D2(*jet) = d2;
294
295 }
296
297 return StatusCode::SUCCESS;
298}
#define M_PI
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
constexpr int pow(int base, int exp) noexcept
WDHK m_dec_ECF1
Energy correlation factors decorators.
SG::ReadHandleKey< xAOD::JetContainer > m_jetKey
WDHK m_dec_pT
Reclustered jets information decorators.
WDHK m_dec_Split12
KtSplittingScale decorators.
std::optional< fastjet::contrib::SoftDrop > m_softdropper
std::optional< fastjet::Filter > m_trimmer
Gaudi::Property< std::vector< std::string > > m_ghostNames
virtual double e() const
The total energy of the particle.
double entries
Definition listroot.cxx:49
time(flags, cells_name, *args, **kw)
float eTot(const U &p)
Accessor utility function for getting the value of Energy.
@ CaloCluster
The object is a calorimeter cluster.
Definition ObjectType.h:39
Jet_v1 Jet
Definition of the current "jet version".
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.

◆ initialize()

StatusCode DerivationFramework::RCJetSubstructureAug::initialize ( )
finaloverridevirtual

Definition at line 29 of file RCJetSubstructureAug.cxx.

30{
31 ATH_MSG_VERBOSE("initialize() ...");
32 ATH_CHECK( m_jetKey.initialize() );
33
34 ATH_CHECK( m_dec_Qw.initialize() );
35 ATH_CHECK( m_dec_Tau1.initialize() );
36 ATH_CHECK( m_dec_Tau2.initialize() );
37 ATH_CHECK( m_dec_Tau3.initialize() );
38 ATH_CHECK( m_dec_Tau4.initialize() );
39 ATH_CHECK( m_dec_Tau21.initialize() );
40 ATH_CHECK( m_dec_Tau32.initialize() );
41 ATH_CHECK( m_dec_Split12.initialize() );
42 ATH_CHECK( m_dec_Split23.initialize() );
43 ATH_CHECK( m_dec_Split34.initialize() );
44 ATH_CHECK( m_dec_ECF1.initialize() );
45 ATH_CHECK( m_dec_ECF2.initialize() );
46 ATH_CHECK( m_dec_ECF3.initialize() );
47 ATH_CHECK( m_dec_ECF4.initialize() );
48 ATH_CHECK( m_dec_C2.initialize() );
49 ATH_CHECK( m_dec_D2.initialize() );
50 ATH_CHECK( m_dec_pT.initialize() );
51 ATH_CHECK( m_dec_m.initialize() );
52 ATH_CHECK( m_dec_NConstits.initialize() );
53 ATH_CHECK( m_dec_eta.initialize() );
54 ATH_CHECK( m_dec_phi.initialize() );
55 ATH_CHECK( m_dec_timing.initialize() );
56
57 // Set up the text-parsing machinery for selectiong the jet directly according to user cuts
58 if (!m_selectionString.empty()) {
59 ATH_CHECK( initializeParser( m_selectionString ));
60 }
61
62 // Define the trimmer
63 m_trimmer.emplace(fastjet::JetDefinition(fastjet::kt_algorithm, m_rclus), fastjet::SelectorPtFractionMin(m_ptfrac));
65
66 return StatusCode::SUCCESS;
67}
#define ATH_MSG_VERBOSE(x)

Member Data Documentation

◆ m_beta

Gaudi::Property<float> DerivationFramework::RCJetSubstructureAug::m_beta {this, "BetaSoft", 1. , "How much to consider angular dependence"}
private

Definition at line 51 of file RCJetSubstructureAug.h.

52{this, "BetaSoft", 1. , "How much to consider angular dependence"};

◆ m_dec_C2

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_C2 { this, "dec_C2", m_jetKey, "C2_" }
private

Definition at line 84 of file RCJetSubstructureAug.h.

84{ this, "dec_C2", m_jetKey, "C2_" };

◆ m_dec_D2

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_D2 { this, "dec_D2", m_jetKey, "D2_" }
private

Definition at line 85 of file RCJetSubstructureAug.h.

85{ this, "dec_D2", m_jetKey, "D2_" };

◆ m_dec_ECF1

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_ECF1 { this, "dec_ECF1", m_jetKey, "ECF1_" }
private

Energy correlation factors decorators.

Definition at line 80 of file RCJetSubstructureAug.h.

80{ this, "dec_ECF1", m_jetKey, "ECF1_" };

◆ m_dec_ECF2

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_ECF2 { this, "dec_ECF2", m_jetKey, "ECF2_" }
private

Definition at line 81 of file RCJetSubstructureAug.h.

81{ this, "dec_ECF2", m_jetKey, "ECF2_" };

◆ m_dec_ECF3

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_ECF3 { this, "dec_ECF3", m_jetKey, "ECF3_" }
private

Definition at line 82 of file RCJetSubstructureAug.h.

82{ this, "dec_ECF3", m_jetKey, "ECF3_" };

◆ m_dec_ECF4

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_ECF4 { this, "dec_ECF4", m_jetKey, "ECF4_" }
private

Definition at line 83 of file RCJetSubstructureAug.h.

83{ this, "dec_ECF4", m_jetKey, "ECF4_" };

◆ m_dec_eta

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_eta { this, "dec_eta", m_jetKey, "eta_" }
private

Definition at line 91 of file RCJetSubstructureAug.h.

91{ this, "dec_eta", m_jetKey, "eta_" };

◆ m_dec_m

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_m { this, "dec_m", m_jetKey, "m_" }
private

Definition at line 89 of file RCJetSubstructureAug.h.

89{ this, "dec_m", m_jetKey, "m_" };

◆ m_dec_NConstits

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_NConstits { this, "dec_NConstits", m_jetKey, "NConstits_" }
private

Definition at line 90 of file RCJetSubstructureAug.h.

90{ this, "dec_NConstits", m_jetKey, "NConstits_" };

◆ m_dec_phi

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_phi { this, "dec_phi", m_jetKey, "phi_" }
private

Definition at line 92 of file RCJetSubstructureAug.h.

92{ this, "dec_phi", m_jetKey, "phi_" };

◆ m_dec_pT

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_pT { this, "dec_pT", m_jetKey, "pT_" }
private

Reclustered jets information decorators.

Definition at line 88 of file RCJetSubstructureAug.h.

88{ this, "dec_pT", m_jetKey, "pT_" };

◆ m_dec_Qw

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_Qw { this, "dec_Qw", m_jetKey, "Qw_" }
private

Qw decorator.

Definition at line 64 of file RCJetSubstructureAug.h.

64{ this, "dec_Qw", m_jetKey, "Qw_" };

◆ m_dec_Split12

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_Split12 { this, "dec_Split12", m_jetKey, "Split12_" }
private

KtSplittingScale decorators.

Definition at line 75 of file RCJetSubstructureAug.h.

75{ this, "dec_Split12", m_jetKey, "Split12_" };

◆ m_dec_Split23

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_Split23 { this, "dec_Split23", m_jetKey, "Split23_" }
private

Definition at line 76 of file RCJetSubstructureAug.h.

76{ this, "dec_Split23", m_jetKey, "Split23_" };

◆ m_dec_Split34

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_Split34 { this, "dec_Split34", m_jetKey, "Split34_" }
private

Definition at line 77 of file RCJetSubstructureAug.h.

77{ this, "dec_Split34", m_jetKey, "Split34_" };

◆ m_dec_Tau1

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_Tau1 { this, "dec_Tau1", m_jetKey, "Tau1_" }
private

Nsubjetiness decorators.

Definition at line 67 of file RCJetSubstructureAug.h.

67{ this, "dec_Tau1", m_jetKey, "Tau1_" };

◆ m_dec_Tau2

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_Tau2 { this, "dec_Tau2", m_jetKey, "Tau2_" }
private

Definition at line 68 of file RCJetSubstructureAug.h.

68{ this, "dec_Tau2", m_jetKey, "Tau2_" };

◆ m_dec_Tau21

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_Tau21 { this, "dec_Tau21", m_jetKey, "Tau21_" }
private

Definition at line 71 of file RCJetSubstructureAug.h.

71{ this, "dec_Tau21", m_jetKey, "Tau21_" };

◆ m_dec_Tau3

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_Tau3 { this, "dec_Tau3", m_jetKey, "Tau3_" }
private

Definition at line 69 of file RCJetSubstructureAug.h.

69{ this, "dec_Tau3", m_jetKey, "Tau3_" };

◆ m_dec_Tau32

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_Tau32 { this, "dec_Tau32", m_jetKey, "Tau32_" }
private

Definition at line 72 of file RCJetSubstructureAug.h.

72{ this, "dec_Tau32", m_jetKey, "Tau32_" };

◆ m_dec_Tau4

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_Tau4 { this, "dec_Tau4", m_jetKey, "Tau4_" }
private

Definition at line 70 of file RCJetSubstructureAug.h.

70{ this, "dec_Tau4", m_jetKey, "Tau4_" };

◆ m_dec_timing

WDHK DerivationFramework::RCJetSubstructureAug::m_dec_timing { this, "dec_timing", m_jetKey, "timing_" }
private

Timing information.

Definition at line 95 of file RCJetSubstructureAug.h.

95{ this, "dec_timing", m_jetKey, "timing_" };

◆ m_ghostNames

Gaudi::Property< std::vector<std::string> > DerivationFramework::RCJetSubstructureAug::m_ghostNames { this, "GhostConstitNames", {"GhostLCTopo"}, "Names of the ghost constituents for substructure computation"}
private

Definition at line 40 of file RCJetSubstructureAug.h.

41{ this, "GhostConstitNames", {"GhostLCTopo"}, "Names of the ghost constituents for substructure computation"};

◆ m_grooming

StringProperty DerivationFramework::RCJetSubstructureAug::m_grooming { this, "Grooming", "", "Name of gromming technic to apply (Trimming or SofDrop)"}
private

Definition at line 44 of file RCJetSubstructureAug.h.

45{ this, "Grooming", "", "Name of gromming technic to apply (Trimming or SofDrop)"};

◆ m_jetKey

SG::ReadHandleKey< xAOD::JetContainer > DerivationFramework::RCJetSubstructureAug::m_jetKey { this, "JetContainerKey", ""}
private

Definition at line 58 of file RCJetSubstructureAug.h.

59{ this, "JetContainerKey", ""};

◆ m_ptfrac

Gaudi::Property<float> DerivationFramework::RCJetSubstructureAug::m_ptfrac {this, "PtFracTrim", 0.03, "pT fraction for retaining subjets"}
private

Definition at line 49 of file RCJetSubstructureAug.h.

50{this, "PtFracTrim", 0.03, "pT fraction for retaining subjets"};

◆ m_R0

Gaudi::Property<float> DerivationFramework::RCJetSubstructureAug::m_R0 {this, "R0Soft", 1. , "Normalization of angular distance, usually the characteristic jet radius (default R0 = 1)"}
private

Definition at line 55 of file RCJetSubstructureAug.h.

56{this, "R0Soft", 1. , "Normalization of angular distance, usually the characteristic jet radius (default R0 = 1)"};

◆ m_rclus

Gaudi::Property<float> DerivationFramework::RCJetSubstructureAug::m_rclus {this, "RClusTrim", 0.3 , "R for reclustering (0 for none)"}
private

Definition at line 47 of file RCJetSubstructureAug.h.

48{this, "RClusTrim", 0.3 , "R for reclustering (0 for none)"};

◆ m_selectionString

StringProperty DerivationFramework::RCJetSubstructureAug::m_selectionString { this, "SelectionString", "", "Selection to apply to jet"}
private

Definition at line 42 of file RCJetSubstructureAug.h.

43{ this, "SelectionString", "", "Selection to apply to jet"};

◆ m_softdropper

std::optional<fastjet::contrib::SoftDrop> DerivationFramework::RCJetSubstructureAug::m_softdropper
private

Definition at line 99 of file RCJetSubstructureAug.h.

◆ m_streamName

StringProperty DerivationFramework::RCJetSubstructureAug::m_streamName { this, "StreamName", "", "Name of the stream" }
private

Definition at line 38 of file RCJetSubstructureAug.h.

39{ this, "StreamName", "", "Name of the stream" };

◆ m_trimmer

std::optional<fastjet::Filter> DerivationFramework::RCJetSubstructureAug::m_trimmer
private

Definition at line 98 of file RCJetSubstructureAug.h.

◆ m_zcut

Gaudi::Property<float> DerivationFramework::RCJetSubstructureAug::m_zcut {this, "ZcutSoft", 1. , "pT fraction for retaining subjets"}
private

Definition at line 53 of file RCJetSubstructureAug.h.

54{this, "ZcutSoft", 1. , "pT fraction for retaining subjets"};

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