ATLAS Offline Software
Loading...
Searching...
No Matches
JetCaloQualityUtils.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7#include "xAODJet/Jet.h"
8
11
12#include "CaloGeoHelpers/CaloSampling.h"
13
14
15using xAOD::Jet;
17
18
19
20namespace {
21
22const int em_calosample[] = { CaloSampling::PreSamplerB, CaloSampling::EMB1, CaloSampling::EMB2, CaloSampling::EMB3,
23 CaloSampling::PreSamplerE, CaloSampling::EME1, CaloSampling::EME2, CaloSampling::EME3,
24 CaloSampling::FCAL0};
25const int had_calosample[] = { CaloSampling::HEC0, CaloSampling::HEC1, CaloSampling::HEC2, CaloSampling::HEC3,
26 CaloSampling::TileBar0, CaloSampling::TileBar1, CaloSampling::TileBar2,
27 CaloSampling::TileGap1, CaloSampling::TileGap2, CaloSampling::TileGap3,
28 CaloSampling::TileExt0, CaloSampling::TileExt1, CaloSampling::TileExt2,
29 CaloSampling::FCAL1, CaloSampling::FCAL2 };
30
31//const double GeV = 1000.0;
32}
33
34namespace jet {
35
36
37 int JetCaloQualityUtils::compute_nLeading(std::vector<double> &cell_energies, const float& e, const float& frac)
38 {
39 std::sort(cell_energies.rbegin(),cell_energies.rend());
40 int counter =0;
41 float sum = 0;
42 for(unsigned int i=0;i<cell_energies.size();i++)
43 {
44 sum += cell_energies[i];
45 counter++;
46 if(sum>frac*e) break;
47 }
48 return counter;
49 }
50
51
52
55
56 double JetCaloQualityUtils::emFraction(const std::vector<float>& e_sampling){
57 double e_EM = 0;
58 for(int i=0; i<9; i++) e_EM += e_sampling[ ::em_calosample[i] ];
59
60 double e_HAD = 0;
61 for(int i=0; i<15; i++) e_HAD += e_sampling[ ::had_calosample[i] ];
62 if( (e_EM==0) || ((e_EM+e_HAD)==0) ) return 0.;
63 return (e_EM / (e_EM+e_HAD));
64 }
65
66
68 {
69
70 const std::vector<float> & einsampling = eperSamplAcc(*jet);
71
72 double e_hec =einsampling[CaloSampling::HEC0]
73 +einsampling[CaloSampling::HEC1]
74 +einsampling[CaloSampling::HEC2]
75 +einsampling[CaloSampling::HEC3];
76
77 double e_jet = jet->jetP4(xAOD::JetEMScaleMomentum).E();
78
79 if(e_jet!=0) return e_hec/e_jet;
80 else return -999;
81 }
82
84 {
85 const std::vector<float>& einsampling = eperSamplAcc(*jet);
86
87 double e_pres = einsampling[CaloSampling::PreSamplerB] + einsampling[CaloSampling::PreSamplerE];
88
89 double e_jet = jet->jetP4(xAOD::JetEMScaleMomentum).E();
90
91 if(e_jet!=0) return e_pres/e_jet;
92 else return -999;
93 }
94
96 {
97 const std::vector<float> & einsampling = eperSamplAcc(*jet);
98
99 double e_tileGap3 =einsampling[CaloSampling::TileGap3];
100 double e_jet = jet->jetP4(xAOD::JetEMScaleMomentum).E();
101
102 if(e_jet!=0) return e_tileGap3/e_jet;
103 else return -999;
104 }
105
106
107 double JetCaloQualityUtils::fracSamplingMax(const Jet* jet, int& SamplingMax)
108 {
109 const std::vector<float> & einsampling = eperSamplAcc(*jet);
110
111 double fracSamplingMax=-999999999.;
112 double sumE_samplings=0.;
113 for ( unsigned int i(0); i < einsampling.size(); ++i )
114 {
115 double e = einsampling[i];
116 sumE_samplings+=e;
117 if (e>fracSamplingMax)
118 {
120 SamplingMax=i;
121 }
122 }
123
124 if(sumE_samplings!=0)
125 fracSamplingMax/=sumE_samplings;
126 else fracSamplingMax=0;
127
128 return fracSamplingMax;
129 }
130
131
132
133
134
135
136
138 {
140 return negE(jet);
141 }
142
143
144 bool JetCaloQualityUtils::isUgly(const Jet* jet,const bool /*recalculateQuantities*/){
145 double fcor = jet->getAttribute<float>(xAOD::JetAttribute::BchCorrCell);
146 double tileGap3f = JetCaloQualityUtils::tileGap3F(jet);
147
148 if (fcor>0.5) return true;
149 if (tileGap3f >0.5) return true;
150
151 return false;
152 }
153
154
155
156
157 // ****************************************************************
158 // JetCalcnLeadingCells *************************************************
159 // ****************************************************************
160
161 bool JetCalcnLeadingCells::setupJet(const Jet* j){
162 m_sumE_cells=0;
163 m_cell_energies.clear();
164 m_cell_energies.reserve(j->numConstituents());
165 return true;
166 }
167
168
170
171 double e = iter->e();
172 if(iter->type() == xAOD::Type::ParticleFlow){
173 e = iter->rawConstituent()->e();
174 }
175 m_cell_energies.push_back(e);
177 return true;
178 }
179
181 std::vector<double> nc_cell_energies(m_cell_energies);
183 }
184
185
186 // ****************************************************************
187 // JetCalcOutOfTimeEnergyFraction *************************************************
188 // ****************************************************************
190 m_sumE=0;
191 m_sumE_OOT=0;
192 return true;
193 }
194
195
197 if( m_sumE== 0.) return -1;
198 return m_sumE_OOT/m_sumE ;
199 }
200
201 bool JetCalcOutOfTimeEnergyFraction::processConstituent(xAOD::JetConstituentVector::iterator& iter){
202
203 //double sum_all(0), sum_time(0);
204
205 double aClusterE = iter->e();
206 if(iter->type() == xAOD::Type::ParticleFlow){
207 aClusterE = iter->rawConstituent()->e();
208 }
209 double aClusterTime = m_constitExtractor->time(iter);//;static_cast<const xAOD::CaloCluster*>(iter->rawConstituent())->time();
210
211
212 if(onlyPosEnergy && aClusterE<0) return true;
213
214 m_sumE += aClusterE;
215 if(fabs(aClusterTime) > timecut) m_sumE_OOT += aClusterE;
216
217 return true;
218 }
219
220 // ****************************************************************
221 // JetCalcTimeCells *************************************************
222 // ****************************************************************
223 bool JetCalcTimeCells::setupJet(const Jet* ){
224 m_time = 0; m_norm = 0;
225 return true;
226 }
227
228
229
230 bool JetCalcTimeCells::processConstituent(xAOD::JetConstituentVector::iterator& iter){
231
232 double thisNorm = iter->e()* iter->e();
233 if(iter->type() == xAOD::Type::ParticleFlow){
234 thisNorm = iter->rawConstituent()->e()*iter->rawConstituent()->e();
235 }
236 m_time += m_constitExtractor->time(iter) * thisNorm; //theClus->time() * thisNorm ;
237 m_norm += thisNorm;
238
239 return true;
240 }
241
242 double JetCalcTimeCells::jetCalculation() const {
243 if (m_norm==0) return 0;
244 return m_time/m_norm ;
245 }
246
247
248 // ****************************************************************
249 // JetCalcAverageLArQualityF *************************************************
250 // ****************************************************************
252 m_qf = 0; m_norm = 0; return true;
253 }
254
255
257 if(m_norm==0) return 0;
258 return m_qf/m_norm;
259 }
260
261 bool JetCalcAverageLArQualityF::processConstituent(xAOD::JetConstituentVector::iterator& iter){
262
263 double e2 = iter->e();
264 if(iter->type() == xAOD::Type::ParticleFlow){
265 e2 = iter->rawConstituent()->e();
266 }
267 e2 = e2*e2;
268
269 m_norm+= e2;
270 double qual(0);
271 if(m_useTile)
273 else
275 m_qf += qual*e2;
276 return true;
277 }
278
279 // ****************************************************************
280 // JetCalcQuality *************************************************
281 // ****************************************************************
282 bool JetCalcQuality::setupJet(const Jet* ){
283 m_totE=0; m_badE=0; return true;
284 }
285
286
287 double JetCalcQuality::jetCalculation() const {
288 if (m_totE==0) return 0;
289 return m_badE/m_totE ;
290 }
291
292
293 bool JetCalcQuality::processConstituent(xAOD::JetConstituentVector::iterator& iter){
294
296
297 if(iter->type() == xAOD::Type::ParticleFlow){
298 m_totE += iter->rawConstituent()->e(); // using iter since it is set at the expected scale by the JetCaloCalculations instance
299 m_badE += f * iter->rawConstituent()->e();
300 }
301 else{
302 m_totE += iter->e();
303 m_badE += f * iter->e();
304 }
305
306 return true;
307 }
308
309
310 // ****************************************************************
311 // JetCalcQualityHEC *************************************************
312 // ****************************************************************
313
314
315
317 double clustHEC = m_constitExtractor->energyHEC(iter);
318 m_totE += clustHEC ;
319
320 double f = m_constitExtractor->moment(iter, xAOD::CaloCluster::BADLARQ_FRAC);
321
322 m_badE += f *clustHEC;
323 return true;
324 }
325
326 // ****************************************************************
327 // JetCalcNegativeEnergy *************************************************
328 // ****************************************************************
330 m_totE =0;
331 return true;
332 }
333
334
336
337 return m_totE;
338
339 }
340
341 bool JetCalcNegativeEnergy::processConstituent(xAOD::JetConstituentVector::iterator& iter){
342
343 double e = iter->e() ; // using iter since it is set at the expected scale by the JetCaloCalculations instance
344 if(iter->type() == xAOD::Type::ParticleFlow){
345 e = iter->rawConstituent()->e();
346 }
347
349 m_totE += (e - epos );
350
351 return true;
352 }
353
354 // ****************************************************************
355 // JetCalcCentroid *************************************************
356 // ****************************************************************
357 bool JetCalcCentroid::setupJet(const Jet* /*j*/){
358 m_totE =0;
360 return true;
361 }
362
363
364 double JetCalcCentroid::jetCalculation() const {
365
366 if (m_totE == 0) return 0;
367 double c_x = m_centroid_x/ m_totE;
368 double c_z = m_centroid_z/ m_totE;
369 double c_y = m_centroid_y/ m_totE;
370
371 return sqrt(c_x*c_x + c_y*c_y+ c_z*c_z);
372 }
373
374 bool JetCalcCentroid::processConstituent(xAOD::JetConstituentVector::iterator& iter){
375
376 double e = iter->e() ; // using iter since it is set at the expected scale by the JetCaloCalculations instance
377 if(iter->type() == xAOD::Type::ParticleFlow){
378 e = iter->rawConstituent()->e();
379 }
380 if(iter->type() == xAOD::Type::FlowElement){
381 e = iter->rawConstituent()->e();
382 }
383
384 m_totE +=e;
385 double x,y,z;
386
390
391 m_centroid_x += x* e;
392 m_centroid_y += y* e;
393 m_centroid_z += z* e;
394
395 return true;
396 }
397
398 // ****************************************************************
399 // JetCalcBadCellsFrac *************************************************
400 // ****************************************************************
401 bool JetCalcBadCellsFrac::setupJet(const Jet* j){
402 m_totE = j->jetP4(xAOD::JetEMScaleMomentum).E() ;
403 m_badE = 0;
404 return true;
405 }
406
407
408 bool JetCalcBadCellsFrac::processConstituent(xAOD::JetConstituentVector::iterator& iter){
410 return true;
411 }
412
414 if (m_totE == 0) return 0;
415 return m_badE / m_totE ;
416 }
417
418
419}
#define y
#define x
#define z
bool processConstituent(xAOD::JetConstituentVector::iterator &iter)
Perform 1 calculation step using 1 constituent.
std::vector< double > m_cell_energies
virtual double jetCalculation() const
return the result of the calculation
virtual bool setupJet(const xAOD::Jet *)=0
virtual bool processConstituent(xAOD::JetConstituentVector::iterator &)
Perform 1 calculation step using 1 constituent.
const CaloConstitHelpers::CaloConstitExtractor * m_constitExtractor
static int compute_nLeading(std::vector< double > &cell_energies, const float &e, const float &frac)
static double emFraction(const std::vector< float > &ePerSampling)
static double hecF(const Jet *jet)
static double fracSamplingMax(const Jet *jet, int &SamplingMax)
static double jetNegativeEnergy(const Jet *jet)
static double presamplerFraction(const Jet *jet)
static double tileGap3F(const Jet *jet)
static bool isUgly(const Jet *jet, const bool recalculateQuantities=false)
@ AVG_TILE_Q
Sum(E_cell_Tile^2 Q_cell_Tile)/Sum(E_cell_Tile^2)
@ AVG_LAR_Q
Sum(E_cell_LAr^2 Q_cell_LAr)/Sum(E_cell_LAr^2)
@ CENTER_Z
Cluster Centroid ( )
@ ENG_BAD_CELLS
Total em-scale energy of bad cells in this cluster.
@ ENG_POS
Total positive Energy of this cluster.
@ BADLARQ_FRAC
Energy fraction of LAr cells with quality larger than a given cut.
@ CENTER_X
Cluster Centroid ( )
@ CENTER_Y
Cluster Centroid ( )
size_t numConstituents() const
Number of constituents in this jets (this is valid even when reading a file where the constituents ha...
Definition Jet_v1.cxx:153
double e2(const xAOD::CaloCluster &cluster)
return the uncorrected cluster energy in 2nd sampling
static const xAOD::JetAttributeAccessor::AccessorWrapper< std::vector< float > > & eperSamplAcc
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
@ ParticleFlow
The object is a particle-flow object.
Definition ObjectType.h:41
@ FlowElement
The object is a track-calo-cluster.
Definition ObjectType.h:52
const AccessorWrapper< T > * accessor(xAOD::JetAttribute::AttributeID id)
Returns an attribute accessor corresponding to an AttributeID.
Jet_v1 Jet
Definition of the current "jet version".
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
@ JetEMScaleMomentum
Definition JetTypes.h:28