ATLAS Offline Software
Loading...
Searching...
No Matches
GetLCOutOfCluster.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5//-----------------------------------------------------------------------
6// File and Version Information:
7// $Id: GetLCOutOfCluster.cxx,v 1.2 2008-11-04 14:33:37 menke Exp $
8//
9// Description: see GetLCOutOfCluster.h
10//
11// Environment:
12// Software developed for the ATLAS Detector at CERN LHC
13//
14// Author List:
15// Sven Menke
16//
17//-----------------------------------------------------------------------
18
19//-----------------------
20// This Class's Header --
21//-----------------------
23
24//---------------
25// C++ Headers --
26//---------------
31
33
34#include "TFile.h"
35#include "TProfile2D.h"
36#include "TString.h"
37#include <iterator>
38#include <cmath>
39
40//###############################################################################
41GetLCOutOfCluster::GetLCOutOfCluster(const std::string& name,
42 ISvcLocator* pSvcLocator)
43 : AthAlgorithm(name, pSvcLocator),
44 m_outputFile(nullptr),
45 m_clusterCollName("CaloTopoClusters"),
50{
51
52 std::vector<Gaudi::Histo1DDef> dims(6);
53 dims[0] = Gaudi::Histo1DDef("side",-1.5,1.5,1);
54 dims[1] = Gaudi::Histo1DDef("|eta|",0.,5.,50);
55 dims[2] = Gaudi::Histo1DDef("phi",-M_PI,M_PI,1);
56 dims[3] = Gaudi::Histo1DDef("log10(E_clus (MeV))",log10(200),log10(1e6),14);
57 dims[4] = Gaudi::Histo1DDef("log10(lambda_clus (mm))",0.0,4.0,14);
58 dims[5] = Gaudi::Histo1DDef("weight",0.,3.0,1);
59
60 mapinsert(dims);
61
62 // Dimensions to use for classification
63 declareProperty("OutOfClusterDimensions",m_dimensionsmap);
64
65 // Name of output file to save histograms in
66 declareProperty("OutputFileName",m_outputFileName);
67
68 // Name of ClusterContainer to use
69 declareProperty("ClusterCollectionName",m_clusterCollName);
70
71 // Name of Samplings to ignore
72 m_invalidSamplingNames.resize(3);
73
74 m_invalidSamplingNames[0] = "PreSamplerB";
75 m_invalidSamplingNames[1] = "PreSamplerE";
76 m_invalidSamplingNames[2] = "TileGap3";
77
78 declareProperty("InvalidSamplings",m_invalidSamplingNames);
79
80 // Normalization type "Const", "Lin", "Log", "NClus"
81 declareProperty("NormalizationType",m_NormalizationType);
82
83 // Classification type "None" for single pion MC or
84 // "ParticleID_EM" for ParticleID based em-type clusters
85 // "ParticleID_HAD" for ParticleID based had-type clusters
86 declareProperty("ClassificationType", m_ClassificationType);
87}
88
89//###############################################################################
90
93
94//###############################################################################
95
97
98 m_outputFile = std::make_unique<TFile>(m_outputFileName.c_str(),"RECREATE");
99 m_outputFile->cd();
100 m_ooc.resize(0);
101
102 mapparse();
103
104 if ( m_NormalizationType == "Lin" ) {
105 ATH_MSG_INFO( "Using weighting proportional to E_calib" );
107 }
108 else if ( m_NormalizationType == "Log" ) {
109 ATH_MSG_INFO( "Using weighting proportional to log(E_calib)" );
111 }
112 else if ( m_NormalizationType == "NClus" ) {
113 ATH_MSG_INFO( "Using weighting proportional to 1/N_Clus_E_calib>0" );
115 }
116 else {
117 ATH_MSG_INFO( "Using constant weighting" );
119 }
120
121 if ( m_ClassificationType == "None" ) {
122 ATH_MSG_INFO( "Expecting single particle input" );
124 }
125 else if ( m_ClassificationType == "ParticleID_EM" ) {
126 ATH_MSG_INFO( "Expecting ParticleID simulation as input -- use EM type clusters only" );
128 }
129 else if ( m_ClassificationType == "ParticleID_HAD" ) {
130 ATH_MSG_INFO( "Expecting ParticleID simulation as input -- use HAD type clusters only" );
132 }
133 else {
134 ATH_MSG_WARNING( " unknown classification type " << m_ClassificationType << " given! Using None instead" );
136 }
137
138 int iside(-1);
139 int ieta(-1);
140 int iphi(-1);
141 int ilogE(-1);
142 int iloglambda(-1);
143 int iweight(-1);
144
145 m_isampmap.resize(3,-1);
146 for(unsigned int idim=0;idim<m_dimensions.size();idim++) {
147 if ( m_dimensions[idim].title() == "side" ) {
148 iside = idim;
149 m_isampmap[0] = iside;
150 }
151 else if ( m_dimensions[idim].title() == "|eta|" )
152 ieta = idim;
153 else if ( m_dimensions[idim].title() == "phi" ) {
154 iphi = idim;
155 m_isampmap[1] = iphi;
156 }
157 else if ( m_dimensions[idim].title() == "log10(E_clus (MeV))" ) {
158 ilogE = idim;
159 m_isampmap[2] = ilogE;
160 }
161 else if ( m_dimensions[idim].title() == "log10(lambda_clus (mm))" )
162 iloglambda = idim;
163 else if ( m_dimensions[idim].title() == "weight" )
164 iweight = idim;
165 }
166 if ( ilogE < 0 || ieta < 0 || iloglambda < 0 || iweight < 0 || iside < 0 ) {
167 ATH_MSG_FATAL(" Mandatory dimension log10E, |eta|, log10lambda or weight missing ...");
168 return StatusCode::FAILURE;
169 }
170 int nside = m_dimensions[iside].bins();
171 int nphi=1;
172 if (iphi>=0) nphi = m_dimensions[iphi].bins();
173 int nlogE = m_dimensions[ilogE].bins();
174 m_ooc.resize(nside*nphi*nlogE,nullptr);
175 for ( int jside=0;jside<nside;jside++) {
176 for ( int jphi=0;jphi<nphi;jphi++) {
177 for(int jlogE=0;jlogE<nlogE;jlogE++) {
178 TString oname("ooc");
179 oname += "_iside_";
180 oname += jside;
181 oname += "_[";
182 oname += m_dimensions[iside].lowEdge();
183 oname += ",";
184 oname += m_dimensions[iside].highEdge();
185 oname += ",";
186 oname += nside;
187 oname += "]";
188 oname += "_iphi_";
189 oname += jphi;
190 oname += "_[";
191 oname += (iphi>=0?m_dimensions[iphi].lowEdge():-1);
192 oname += ",";
193 oname += (iphi>=0?m_dimensions[iphi].highEdge():-1);
194 oname += ",";
195 oname += nphi;
196 oname += "]";
197 oname += "_ilogE_";
198 oname += jlogE;
199 oname += "_[";
200 oname += m_dimensions[ilogE].lowEdge();
201 oname += ",";
202 oname += m_dimensions[ilogE].highEdge();
203 oname += ",";
204 oname += nlogE;
205 oname += "]";
206 int iO = jlogE*nphi*nside+jphi*nside+jside;
207 m_ooc[iO] = new TProfile2D(oname,oname,
208 m_dimensions[ieta].bins(),
209 m_dimensions[ieta].lowEdge(),
210 m_dimensions[ieta].highEdge(),
211 m_dimensions[iloglambda].bins(),
212 m_dimensions[iloglambda].lowEdge(),
213 m_dimensions[iloglambda].highEdge(),
214 m_dimensions[iweight].lowEdge(),
215 m_dimensions[iweight].highEdge());
216 m_ooc[iO]->SetXTitle("|#eta_{clus}|");
217 m_ooc[iO]->SetYTitle("log_{10}(#lambda_{clus} (mm))");
218 m_ooc[iO]->SetZTitle("E_{out of cluster} / E_{clus}^{EM-no-PS/Gap3} / Isolation");
219 }
220 }
221 }
222 //--- check sampling names to use exclude in correction
223 for (const std::string& name : m_invalidSamplingNames) {
224 int theSampling(CaloSampling::Unknown);
225 for (unsigned int jsamp = 0;jsamp< CaloSampling::Unknown; jsamp++) {
227 theSampling = jsamp;
228 break;
229 }
230 }
231 if ( theSampling == CaloSampling::Unknown ) {
232 msg(MSG::ERROR) << "Calorimeter sampling "
233 << name
234 << " is not a valid Calorimeter sampling name and will be ignored! "
235 << "Valid names are: ";
236 for (unsigned int jsamp = 0;jsamp< CaloSampling::Unknown; jsamp++) {
238 if ( jsamp < CaloSampling::Unknown-1)
239 msg() << ", ";
240 else
241 msg() << ".";
242 }
243 msg() << endmsg;
244 }
245 else {
246 m_invalidSamplings.insert(theSampling);
247 }
248 }
249
250 msg(MSG::INFO) << "Samplings to exclude from the out-of-cluster weighting:";
251 for (const std::string& name : m_invalidSamplingNames)
252 msg() << " " << name;
253 msg() << endmsg;
254
255 ATH_CHECK( m_clusterCollName.initialize() );
256
257 return StatusCode::SUCCESS;
258}
259
260//###############################################################################
261
263{
264 ATH_MSG_INFO( "Writing out histograms" );
265 m_outputFile->cd();
266 for(unsigned int i=0;i<m_ooc.size();i++) {
267 m_ooc[i]->Write();
268 }
269 m_outputFile->Close();
270
271 return StatusCode::SUCCESS;
272}
273
274//###############################################################################
275
277{
279
280 // total calib hit energy of all clusters
281 double eCalibTot(0.);
282 double nClusECalibGt0 = 0.0;
283
284 for (const xAOD::CaloCluster* theCluster : *cc) {
285 double eC=999;
286 if (!theCluster->retrieveMoment(xAOD::CaloCluster::ENG_CALIB_TOT,eC)) {
287 ATH_MSG_ERROR( "Failed to retrieve cluster moment ENG_CALIB_TOT" );
288 return StatusCode::FAILURE;
289 }
291 double emFrac=-999;
292 if (!theCluster->retrieveMoment(xAOD::CaloCluster::ENG_CALIB_FRAC_EM,emFrac)){
293 ATH_MSG_ERROR( "Failed to retrieve cluster moment ENG_CALIB_FAC_EM" );
294 return StatusCode::FAILURE;
295 }
297 eC = 0;
299 eC = 0;
300 }
301 eCalibTot += eC;
302 if ( eC > 0 ) {
303 nClusECalibGt0++;
304 }
305 }
306
307 if ( eCalibTot > 0 && nClusECalibGt0 > 0 ) {
308 const double inv_eCalibTot = 1. / eCalibTot;
309 const double inv_nClusECalibGt0 = 1. / nClusECalibGt0;
310 for (const xAOD::CaloCluster* pClus : *cc) {
311 double eng = pClus->e();
312
314 double emFrac=-999;
315 if (!pClus->retrieveMoment(xAOD::CaloCluster::ENG_CALIB_FRAC_EM,emFrac)){
316 ATH_MSG_ERROR( "Failed to retrieve cluster moment ENG_CALIB_FAC_EM");
317 return StatusCode::FAILURE;
318 }
320 continue;
322 continue;
323 }
324
325 // subtract the samplings to ignore from eng
326 for (int sampling : m_invalidSamplings) {
327 eng -= pClus->eSample((CaloSampling::CaloSample)(sampling));
328 }
329
330 if ( eng > 0 ) {
331
332 int iside = 0;
333 int iphi = 0;
334 int ilogE = 0;
335 int nside = 1;
336 int nphi = 1;
337 int nlogE = 1;
338
339 if ( m_isampmap[0] >= 0 ) {
340 const Gaudi::Histo1DDef & hd = m_dimensions[m_isampmap[0]];
341 nside = hd.bins();
342 iside = (int)(nside*(((pClus->eta()<0?-1.0:1.0) - hd.lowEdge())
343 /(hd.highEdge()-hd.lowEdge())));
344 if ( iside < 0 || iside > nside-1 ) {
345 ATH_MSG_WARNING( " Side index out of bounds " <<
346 iside << " not in [0," << nside-1 << "]" );
347 iside = -1;
348 }
349 }
350
351 if ( m_isampmap[1] >= 0 ) {
352 const Gaudi::Histo1DDef & hd = m_dimensions[m_isampmap[1]];
353 nphi = hd.bins();
354 iphi = (int)(nphi*((pClus->phi() - hd.lowEdge())
355 /(hd.highEdge()-hd.lowEdge())));
356 if ( iphi < 0 || iphi > nphi-1 ) {
357 ATH_MSG_WARNING( " Phi index out of bounds " <<
358 iphi << " not in [0," << nphi-1 << "]" );
359 iphi = -1;
360 }
361 }
362
363 const Gaudi::Histo1DDef & hd = m_dimensions[m_isampmap[2]];
364 nlogE = hd.bins();
365 ilogE = (int)(nlogE*((log10(eng) - hd.lowEdge())
366 /(hd.highEdge()-hd.lowEdge())));
367 if ( ilogE >= 0 && ilogE < nlogE ) {
368 double lamb,eout,etot,isol;
369 if (!pClus->retrieveMoment(xAOD::CaloCluster::CENTER_LAMBDA,lamb) ||
370 !pClus->retrieveMoment(xAOD::CaloCluster::ENG_CALIB_OUT_L,eout) ||
371 !pClus->retrieveMoment(xAOD::CaloCluster::ENG_CALIB_TOT,etot) ||
372 !pClus->retrieveMoment(xAOD::CaloCluster::ISOLATION,isol)) {
373 ATH_MSG_ERROR( "Failed to retrieve a cluster moment (CENTER_LAMBDA,ENG_CALIB_OUT,ENG_CALIB_TOT,ISOLATION)" );
374 return StatusCode::FAILURE;
375 }
376
377 if ( lamb > 0 &&
378 etot > 0 &&
379 isol > 0.5 )
380 {
381 int iO = ilogE*nphi*nside+iphi*nside+iside;
382 if (m_ooc[iO]) {
383 double norm = 0.0;
385 norm = etot*inv_eCalibTot;
386 }
388 // cluster has to have at least 1% of the calib hit E
389 norm = log10(etot*inv_eCalibTot)+2.0;
390 }
392 norm = inv_nClusECalibGt0;
393 }
394 else {
395 norm = 1.0;
396 }
397 if ( norm > 0 ) {
398 m_ooc[iO]->Fill(fabs(pClus->eta()),log10(lamb),eout/eng/isol,norm);
399 }
400 }
401 }
402 }
403 }
404 }
405 }
406
407 return StatusCode::SUCCESS;
408}
409
410void GetLCOutOfCluster::mapinsert(const std::vector<Gaudi::Histo1DDef> & dims) {
411 for (unsigned int i=0;i<dims.size();i++) {
412 m_dimensionsmap[dims[i].title()] = dims[i];
413 }
414}
415
417
418 for (std::pair<const std::string, Gaudi::Histo1DDef>& p : m_dimensionsmap) {
419 m_dimensions.push_back(p.second);
420 ATH_MSG_DEBUG(" New Dimension: "
421 << p.second.title() << ", [" << p.second.lowEdge()
422 << ", " << p.second.highEdge()
423 << ", " << p.second.bins()
424 << "]");
425 }
426}
427
428//###############################################################################
#define M_PI
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Helpers for checking error return status codes and reporting errors.
static const std::vector< std::string > bins
Handle class for reading from StoreGate.
AthAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
MsgStream & msg() const
static const std::string & getSamplingName(const CaloSampling::CaloSample theSample)
Returns a string (name) for each CaloSampling.
std::string m_ClassificationType
string to choose different classification types
std::vector< std::string > m_invalidSamplingNames
vector of names of the calorimeter samplings not to use when applying the out-of-cluster weights.
std::string m_outputFileName
Name of the output file to save histograms in.
SG::ReadHandleKey< xAOD::CaloClusterContainer > m_clusterCollName
Name of the CaloClusterContainer to use.
std::set< int > m_invalidSamplings
actual set of samplings to be ignored for out-of-cluster weights
std::vector< TProfile2D * > m_ooc
Vector of actual histograms.
std::vector< Gaudi::Histo1DDef > m_dimensions
definition of all dimensions used for out-of-cluster corrections
virtual StatusCode finalize()
void mapinsert(const std::vector< Gaudi::Histo1DDef > &dims)
virtual StatusCode initialize()
std::string m_NormalizationType
string to choose different normalization types
std::vector< int > m_isampmap
Vector of indices in m_dimensions.
GetLCOutOfCluster(const std::string &name, ISvcLocator *pSvcLocator)
virtual StatusCode execute()
std::map< std::string, Gaudi::Histo1DDef > m_dimensionsmap
property to set all dimensions introduced above
std::unique_ptr< TFile > m_outputFile
Output file to save histograms in.
@ ENG_CALIB_OUT_L
Attached Calibration Hit energy outside clusters but inside the calorimeter with loose matching (Angl...
@ CENTER_LAMBDA
Shower depth at Cluster Centroid.
@ ENG_CALIB_FRAC_EM
Calibration Hit energy inside the cluster caused by e/gamma/pi0.
@ ENG_CALIB_TOT
Calibration Hit energy inside the cluster.
@ ISOLATION
Energy weighted fraction of non-clustered perimeter cells.
@ PARTICLEID_EM
Definition GetLCDefs.h:25
@ PARTICLEID_HAD
Definition GetLCDefs.h:25
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.