ATLAS Offline Software
Loading...
Searching...
No Matches
CaloRingsBuilder.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#include "CaloRingsBuilder.h"
7
9
10// Cell includes:
12#include "CaloGeoHelpers/CaloSampling.h"
14
15// RingSet includes:
18
19// CaloRings includes:
22
23// Ringer Conf include:
25
26// Other xAOD includes:
27#include "xAODBase/IParticle.h"
29
30// STL
31#include <algorithm>
32#include <cfloat>
33#include <cmath>
34#include <sstream>
35
36namespace Ringer {
37
38// =================================================================================
40 const std::string& name,
41 const ::IInterface* parent)
42 : ::AthAlgTool(type, name, parent),
43 m_rsCont(nullptr),
44 m_crCont(nullptr),
46{
47 declareInterface<ICaloRingsBuilder>(this);
48}
49
50// =====================================================================================
52
53// =====================================================================================
55{
56 ATH_MSG_DEBUG("Initializing " << name() );
57
58 m_nRingSets = m_nRings.size();
59
60 auto itr = m_layers.begin();
61
62 // Build RingSets configuration:
63 for (size_t rsConfIdx = 0; rsConfIdx < m_nRingSets; ++rsConfIdx) {
64
65 const auto rsNLayers = m_nLayers[rsConfIdx];
66
67 auto end_itr = itr + rsNLayers;
68
69 const auto& caloSampleItr = reinterpret_cast<
70 std::vector<CaloSampling::CaloSample>::iterator&
71 >(itr);
72 const auto& caloSampleEndItr = reinterpret_cast<
73 std::vector<CaloSampling::CaloSample>::iterator&
74 >(end_itr);
75
76 std::vector<CaloSampling::CaloSample> rsLayers( caloSampleItr ,
77 caloSampleEndItr);
78
79 itr += rsNLayers;
80
82 m_nRings[rsConfIdx],
83 rsLayers,
84 m_etaWidth[rsConfIdx], m_phiWidth[rsConfIdx],
88 );
89
90 // Build our raw configuration structure:
91 m_rsRawConfCol.push_back(std::move(rawConf));
92 }
93
94 // We have finished filling the main raw configuration properties, now we add
95 // it bounderies:
96 try {
98 } catch ( const std::runtime_error &e) {
99 ATH_MSG_ERROR("Could not add collection bounderies due to: " << e.what() );
100 ATH_MSG_ERROR("RawConfCollection is: ");
101 std::ostringstream str;
103 ATH_MSG_ERROR(str.str());
104 return StatusCode::FAILURE;
105 }
106
107 // Print our collection
108 if (msgLevel() <= MSG::DEBUG){
109 std::ostringstream str;
111 ATH_MSG_DEBUG(str.str());
112 }
113
114 // This will check that the properties were initialized properly
115 // by job configuration.
116 ATH_CHECK( m_crContName.initialize() );
117 ATH_CHECK( m_rsContName.initialize() );
118 ATH_CHECK( m_cellsContName.initialize() );
119 ATH_CHECK( m_caloMgrKey.initialize() );
120
121 return StatusCode::SUCCESS;
122}
123
124// =====================================================================================
126{
127 return StatusCode::SUCCESS;
128}
129
130// =====================================================================================
132 , xAOD::RingSetContainer* rsCont
133 , std::size_t nReserve )
134{
135 if ( crCont && rsCont ){
136 m_crCont = crCont;
137 m_rsCont = rsCont;
138 } else {
139 ATH_MSG_ERROR( "Attempted to set CaloRingsContainer and/or RingSetContainer to invalid pointers");
140 return StatusCode::FAILURE;
141 }
142 // Reserve container space if required:
143 if (nReserve) {
144 // Reserve one CaloRings per particle
145 m_crCont->reserve( nReserve );
146 // We need to reserve more space for the RingSet container, there will be
147 // the number of RawConfig available in our raw configuration collection.
148 m_rsCont->reserve( nReserve * m_nRingSets );
149 }
150
151 return StatusCode::SUCCESS;
152}
153
154// =====================================================================================
157{
158 double et(0.);
159 const double eta2 = std::fabs(cluster.etaBE(2));
160 const double energy = cluster.e();
161 if ( eta2 < 999.) {
162 const double cosheta = std::cosh(eta2);
163 et = (cosheta != 0.) ? energy /cosheta : 0.;
164 }
165 if ( et > m_minEnergy )
166 {
167 return executeTemp(cluster, clRings);
168 } else {
169 ATH_MSG_DEBUG( "Skipping cluster with low energy (" << et << " MeV) .");
170 return StatusCode::SUCCESS;
171 }
172}
173
174// =====================================================================================
176 const xAOD::IParticle &particle,
178{
179 double et = particle.pt();
180 if ( et > m_minEnergy )
181 {
182 return executeTemp(particle, clRings);
183 } else {
184 ATH_MSG_DEBUG( "Skipping particle with low energy (" << et << " MeV) .");
185 return StatusCode::SUCCESS;
186 }
187}
188
189// =====================================================================================
190// Local execute
191// =====================================================================================
192template<typename T>
194 const T &input,
196{
197
198 ATH_MSG_DEBUG("Entering "<< name() << " execute.");
199
200 // Create structure to hold rings:
201 xAOD::CaloRings *clRings = new xAOD::CaloRings();
202
203 // Add the CaloRings to the container:
204 m_crCont->push_back(clRings);
205
206 // Set elementLink reference to created CaloRings:
207 clRingsLink.toContainedElement(*m_crCont, clRings);
208
209 // If not using shower barycenter, we need to reset last valid seed to avoid
210 // any possible issue:
212 m_lastValidSeed = AtlasGeoPoint(input.eta(),input.phi());
213 }
214
215 // Build this CaloRings RingSets:
216 for ( const auto &rawConf : m_rsRawConfCol ) {
217
218 // Create new RingSet and add it to the container:
219 xAOD::RingSet* rs = new xAOD::RingSet( rawConf.nRings );
220 m_rsCont->push_back(rs);
221
222 // Get RingSet seed:
223 AtlasGeoPoint seed;
224 CHECK( getRingSetSeed( rawConf, input, seed ) );
225
226 // Build it:
227 CHECK( buildRingSet( rawConf, seed, rs ) );
228
229 // Get the ElementLink and add it to our RingSet holder:
231 clRings->addRingSetEL( rsEL );
232 }
233
234 // Print CaloRings with DEBUG level:
235 if (msgLevel() <= MSG::DEBUG){
236 std::ostringstream str;
237 clRings->print(str);
238 ATH_MSG_DEBUG(str.str());
239 }
240
241 return StatusCode::SUCCESS;
242}
243
244// =================================================================================
247 const xAOD::CaloCluster &cluster,
248 AtlasGeoPoint &seed )
249{
251
252 seed.setEta( cluster.eta() );
253 seed.setPhi( cluster.phi() );
254
255 return StatusCode::SUCCESS;
256
257 } else {
258
259 bool foundValid = false, foundMultipleValid = false;
260
261 for ( const auto layer : rawConf.layers ) {
262
263 AtlasGeoPoint seedCandidate(
264 cluster.etaSample( layer ),
265 cluster.phiSample( layer )
266 );
267
268 ATH_MSG_DEBUG( "This layer (" << CaloSampling::getSamplingName(layer) <<
269 ") seedCandidate is: (" <<
270 seedCandidate.eta() << "," <<
271 seedCandidate.phi() << ")" );
272
273 if ( seedCandidate.isValid() ){
274 m_lastValidSeed = seedCandidate;
275 if ( foundValid )
276 foundMultipleValid = true;
277 foundValid = true;
278 }
279 }
280
281 // If we get here, set it to the last valid seed:
282 seed = m_lastValidSeed;
283 if ( foundMultipleValid ){
284 ATH_MSG_WARNING( "Found multiple valid seeds. Set to last valid candidate.");
285 }
286 return StatusCode::SUCCESS;
287 }
288}
289
290// =================================================================================
292 const xAOD::RingSetConf::RawConf &/*rawConf*/,
293 const xAOD::IParticle &part,
294 AtlasGeoPoint &seed)
295{
296
297 seed.setEta( part.eta() );
298 seed.setPhi( part.phi() );
299
300 return StatusCode::SUCCESS;
301}
302
303// =================================================================================
304// StatusCode CaloRingsBuilder::getRingSetSeed(
305// const xAOD::RingSetConf::RawConf &/*rawConf*/,
306// const xAOD::Jet_v1 &jet,
307// AtlasGeoPoint &seed)
308// {
309
310// seed.setEta( jet.eta() );
311// seed.setPhi( jet.phi() );
312
313// return StatusCode::SUCCESS;
314// }
315
316// =================================================================================
319 const AtlasGeoPoint &seed,
321{
322 // Get this RingSet size:
323 const auto nRings = rawConf.nRings;
324
325 // Retrieve CaloCells
327 // check is only used for serial running; remove when MT scheduler used
328 if(!cellsCont.isValid()) {
329 ATH_MSG_FATAL("Failed to retrieve "<< m_cellsContName.key());
330 return StatusCode::FAILURE;
331 }
332
334 const CaloDetDescrManager* caloMgr=*caloMgrHandle;
335
336 CaloCellList cells(caloMgr,cellsCont.ptr() );
337
338 // loop over cells
339 for ( const int layer : rawConf.layers) { // We use int here because the
340 // cells.select() method needs int as param
341 cells.select(seed.eta(), seed.phi(), m_cellMaxDEtaDist, m_cellMaxDPhiDist, layer );
342 for ( const CaloCell *cell : cells ) {
343
344 unsigned int ringNumber(0);
345
346 // calculate the normalised difference in eta
347 const float deltaEta = fabs(
348 (cell->eta() - seed.eta())
349 ) / rawConf.etaWidth;
350 // calculate the normalised difference in phi
351 const float deltaPhi = fabs(
352 CaloPhiRange::diff(cell->phi(), seed.phi())
353 ) / rawConf.phiWidth;
354 // The biggest difference indicates the ring number (we are using
355 // squared-shape rings)
356 const float deltaGreater = std::max(deltaEta, deltaPhi);
357
358 // Round to nearest integer:
359 ringNumber = static_cast<unsigned int>(std::floor(deltaGreater + .5));
360 if ( ringNumber < nRings ){
362 rs->at(ringNumber) += cell->energy()/cosh(cell->eta());
363 }else{
364 rs->at(ringNumber) += cell->energy();
365 }
366 }
367 }
368 }
369
370 return StatusCode::SUCCESS;
371}
372
373
374} // namespace Ringer
Scalar deltaPhi(const MatrixBase< Derived > &vec) const
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Helpers for checking error return status codes and reporting errors.
#define CHECK(...)
Evaluate an expression and check for errors.
float et(const xAOD::jFexSRJetRoI *j)
static Double_t rs
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
float eta() const
bool isValid() const
float phi() const
Data object for each calorimeter readout cell.
Definition CaloCell.h:57
This class provides the client interface for accessing the detector description information common to...
static double diff(double phi1, double phi2)
simple phi1 - phi2 calculation, but result is fixed to respect range.
static std::string getSamplingName(CaloSample theSample)
Returns a string (name) for each CaloSampling.
virtual StatusCode buildRingSet(const xAOD::RingSetConf::RawConf &rawConf, const AtlasGeoPoint &seed, xAOD::RingSet *rs)
main method where the RingSets are built.
CaloRingsBuilder(const std::string &type, const std::string &name, const IInterface *parent)
Default constructor.
~CaloRingsBuilder()
Destructor.
SG::WriteHandleKey< xAOD::RingSetContainer > m_rsContName
Name of RingSetContainer on Event StoreGate.
const xAOD::RingSetConf::RawConfCollection & rawConf() override
Extra methods:
xAOD::RingSetConf::RawConfCollection m_rsRawConfCol
holds each RingSet configuration (filled at initialize)
virtual StatusCode execute(const xAOD::IParticle &particle, ElementLink< xAOD::CaloRingsContainer > &clRingsLink) override
build CaloRings for IParticle
StatusCode getRingSetSeed(const xAOD::RingSetConf::RawConf &conf, const xAOD::CaloCluster &cluster, AtlasGeoPoint &seed)
Fill RingSet seed for CaloCluster.
SG::WriteHandleKey< xAOD::CaloRingsContainer > m_crContName
Fill RingSet seed for IParticle.
Gaudi::Property< float > m_minEnergy
Minimum particle energy to build rings (GeV)
size_t m_nRingSets
hold the number of RingSets we are building for each CaloRings
Gaudi::Property< std::vector< float > > m_etaWidth
Width of the ring in eta.
Gaudi::Property< bool > m_doTransverseEnergy
Switch to use raw cell energy instead ET.
StatusCode executeTemp(const T &input, ElementLink< xAOD::CaloRingsContainer > &crEL)
Tool protected methods:
Gaudi::Property< std::vector< unsigned int > > m_nLayers
Number of calorimeter layers in each ringset.
Gaudi::Property< float > m_cellMaxDEtaDist
Maximum cell distance in eta to seed.
xAOD::RingSetContainer * m_rsCont
Tool props (non configurables):
SG::ReadCondHandleKey< CaloDetDescrManager > m_caloMgrKey
virtual StatusCode preExecute(xAOD::CaloRingsContainer *crCont, xAOD::RingSetContainer *rsCont, const std::size_t nReserve=0) override
method for working on containers
Gaudi::Property< std::vector< float > > m_phiWidth
Width of the ring in phi.
xAOD::CaloRingsContainer * m_crCont
Create and hold CaloRingsContainer for each event.
Gaudi::Property< std::vector< int > > m_layers
Calorimeter layers in each ringset.
SG::ReadHandleKey< CaloCellContainer > m_cellsContName
Name of CaloCellContainer.
AtlasGeoPoint m_lastValidSeed
last valid RingSet seed
Gaudi::Property< std::vector< unsigned int > > m_nRings
Number of rings in a ringset.
Gaudi::Property< bool > m_useShowShapeBarycenter
Switch to use shower barycenter seed for each RingSets.
virtual StatusCode finalize() override
finalize method
virtual StatusCode initialize() override
Tool main methods:
Gaudi::Property< float > m_cellMaxDPhiDist
Maximum cell distance in phi to seed.
const_pointer_type ptr()
Dereference the pointer.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
float phiSample(const CaloSample sampling) const
Retrieve barycenter in a given sample.
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double e() const
The total energy of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
float etaSample(const CaloSample sampling) const
Retrieve barycenter in a given sample.
float etaBE(const unsigned layer) const
Get the eta in one layer of the EM Calo.
void print(std::ostream &stream) const
Print-out methods:
void addRingSetEL(const ElementLink< RingSetContainer_v1 > &rsEL)
Add ElementLink to holden vector.
Class providing the definition of the 4-vector interface.
static Ringer::CalJointLayer whichLayer(const std::vector< CaloSampling::CaloSample > &layers)
static void print(const RawConf &raw, std::ostream &stream)
Prints rawConf.
static Ringer::CalJointSection whichSection(const std::vector< CaloSampling::CaloSample > &layers)
static void addRawConfColBounderies(RawConfCollection &clRingsConf)
Add to RawConfCollection its JointLayer/JointSection bounderies.
Namespace dedicated for Ringer utilities.
RingSet_v1 RingSet
Definition of the current "RingSet version".
Definition RingSet.h:15
CaloRings_v1 CaloRings
Definition of the current "CaloRings version".
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
RingSetContainer_v1 RingSetContainer
Definition of the current "RingSet container version".
CaloRingsContainer_v1 CaloRingsContainer
Definition of the current "CaloRings container version".
Extra patterns decribing particle interation process.