ATLAS Offline Software
Loading...
Searching...
No Matches
CaloClusterMomentsMaker Class Reference

Calculate moments for CaloCluster objects. More...

#include <CaloClusterMomentsMaker.h>

Inheritance diagram for CaloClusterMomentsMaker:
Collaboration diagram for CaloClusterMomentsMaker:

Public Member Functions

 CaloClusterMomentsMaker (const std::string &type, const std::string &name, const IInterface *parent)
virtual StatusCode execute (const EventContext &ctx, xAOD::CaloClusterContainer *theClusColl) const override final
 Execute on an entire collection of clusters.
virtual StatusCode initialize () override
virtual StatusCode finalize () override
virtual StatusCode execute (xAOD::CaloClusterContainer *collection) final
 Execute on an entire collection of clusters.
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysInitialize () override
 Perform system initialization for an algorithm.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const
 DeclareInterfaceID (CaloClusterCollectionProcessor, 1, 0)

Protected Member Functions

void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

std::vector< std::string > m_momentsNames
 vector holding the input list of names of moments to calculate.
std::vector< xAOD::CaloCluster::MomentTypem_validMoments
 set of moments which will be calculated.
const CaloCell_IDm_calo_id
double m_maxAxisAngle
 the maximal allowed deviation from the IP-to-ClusterCenter-axis.
double m_minRLateral
 the minimal \(r\) in the definition of the Lateral moment
double m_minLLongitudinal
 the minimal \(\lambda\) in the definition of the Longitudinal moment
double m_minBadLArQuality
 the minimal cell quality in the LAr for declaring a cell bad
bool m_calculateSignificance
 Set to true if significance moments are need.
bool m_calculateIsolation
 Set to true if cluster isolation is to be calculated.
bool m_calculateLArHVFraction
 Set to true to calculate E and N of cells affected by LAr HV corrections.
bool m_twoGaussianNoise
 if set to true use 2-gaussian noise description for TileCal
ToolHandle< CaloDepthToolm_caloDepthTool
SG::ReadCondHandleKey< CaloDetDescrManagerm_caloMgrKey
SG::ReadCondHandleKey< CaloNoisem_noiseCDOKey {this,"CaloNoiseKey","totalNoise","SG Key of CaloNoise data object"}
 Key of the CaloNoise Conditions data object.
ToolHandle< ILArHVFractionm_larHVFraction
std::string m_momentsNamesAOD
 Not used anymore (with xAOD), but required when configured from COOL.
bool m_absOpt
 if set to true use abs E value of cells to calculate cluster moments
bool m_secondTime = { false }
 Retrieve second moment of cell times and store as moment.
bool m_nCellsPerSampling = { false }
 store number of cells per sampling layer as moment
double m_etaInnerWheel = { 2.52 }
 Transition from outer to inner wheel in EME2.
Gaudi::Property< bool > m_useGPUCriteria {this, "UseGPUCriteria", false, "Adopt a set of criteria that is consistent with the GPU implementation."}
StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Detailed Description

Calculate moments for CaloCluster objects.

Author
Sven Menke menke.nosp@m.@mpp.nosp@m.mu.mp.nosp@m.g.de
Peter Loch loch@.nosp@m.phys.nosp@m.ics.a.nosp@m.rizo.nosp@m.na.ed.nosp@m.u
Date
23-March-2021

This is a CaloClusterCollectionProcessor which can be plugged into a CaloClusterMaker for calculating moments for each CaloCluster. Typically a cluster moment of a certain degree \(n\) over a variable \(x\) defined for each CaloCell member of a CaloCluster is defined as:

\[ \langle x^n\rangle = \frac{1}{E_{\rm norm}} \times \!\!\!\sum\limits_{\{{\rm cell}|E_{\rm cell}>0\}}{\!\!\!E_{\rm cell}\, x^n},\]

with

\[ E_{\rm norm} = \!\!\!\sum\limits_{\{{\rm cell}|E_{\rm cell}>0\}}{\!\!\!E_{\rm cell}}. \]

Note that only cells with positive energy are used in this definition. Common variables to calculate first and second moments of are \(\phi\), \(\eta\), and radial and longitudinal distances from the shower axis and the shower center, respectively.

Since
23-March-2021: second moment of cell time distribution is calculated
Author
Sven Menke menke.nosp@m.@mpp.nosp@m.mu.mp.nosp@m.g.de
Date
28-February-2005

This is a CaloClusterCollectionProcessor which can be plugged into a CaloClusterMaker for calculating moments for each CaloCluster. Typically a cluster moment of a certain degree \(n\) over a variable \(x\) defined for each CaloCell member of a CaloCluster is defined as:

\[ \langle x^n\rangle = \frac{1}{E_{\rm norm}} \times \!\!\!\sum\limits_{\{{\rm cell}|E_{\rm cell}>0\}}{\!\!\!E_{\rm cell}\, x^n},\]

with

\[ E_{\rm norm} = \!\!\!\sum\limits_{\{{\rm cell}|E_{\rm cell}>0\}}{\!\!\!E_{\rm cell}}. \]

Note that only cells with positive energy are used in this definition. Common variables to calculate first and second moments of are \(\phi\), \(\eta\), and radial and longitudinal distances from the shower axis and the shower center, respectively.

Definition at line 49 of file CaloClusterMomentsMaker.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< AlgTool > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ CaloClusterMomentsMaker()

CaloClusterMomentsMaker::CaloClusterMomentsMaker ( const std::string & type,
const std::string & name,
const IInterface * parent )

Definition at line 141 of file CaloClusterMomentsMaker.cxx.

144 : AthAlgTool(type, name, parent),
145 m_calo_id(nullptr),
146 m_maxAxisAngle(20*deg),
147 m_minRLateral(4*cm),
149 m_minBadLArQuality(4000),
153 m_twoGaussianNoise(false),
154 m_caloDepthTool("CaloDepthTool",this),
155 m_larHVFraction("LArHVFraction",this),
156 m_absOpt(false)
157{
158 declareInterface<CaloClusterCollectionProcessor> (this);
159 // Name(s) of Moments to calculate
160 declareProperty("MomentsNames",m_momentsNames);
161
162 // Maximum allowed angle between shower axis and the vector pointing
163 // to the shower center from the IP in degrees. This property is needed
164 // to protect against cases where all significant cells are in one sampling
165 // and the shower axis can thus not be defined.
166 declareProperty("MaxAxisAngle",m_maxAxisAngle);
167 declareProperty("MinRLateral",m_minRLateral);
168 declareProperty("MinLLongitudinal",m_minLLongitudinal);
169 declareProperty("MinBadLArQuality",m_minBadLArQuality);
170 // Use 2-gaussian noise for Tile
171 declareProperty("TwoGaussianNoise",m_twoGaussianNoise);
172 declareProperty("LArHVFraction",m_larHVFraction,"Tool Handle for LArHVFraction");
173 // Not used anymore (with xAOD), but required when configured from COOL.
174 declareProperty("AODMomentsNames",m_momentsNamesAOD);
175 // Use weighting of neg. clusters option?
176 declareProperty("WeightingOfNegClusters", m_absOpt);
177 // Set eta boundary for transition from outer to inner wheel in EME2
178 declareProperty("EMECAbsEtaWheelTransition",m_etaInnerWheel);
179}
#define deg
AthAlgTool()
Default constructor:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
double m_maxAxisAngle
the maximal allowed deviation from the IP-to-ClusterCenter-axis.
bool m_calculateIsolation
Set to true if cluster isolation is to be calculated.
double m_minRLateral
the minimal in the definition of the Lateral moment
bool m_calculateSignificance
Set to true if significance moments are need.
bool m_calculateLArHVFraction
Set to true to calculate E and N of cells affected by LAr HV corrections.
ToolHandle< CaloDepthTool > m_caloDepthTool
std::string m_momentsNamesAOD
Not used anymore (with xAOD), but required when configured from COOL.
bool m_absOpt
if set to true use abs E value of cells to calculate cluster moments
double m_minLLongitudinal
the minimal in the definition of the Longitudinal moment
std::vector< std::string > m_momentsNames
vector holding the input list of names of moments to calculate.
double m_etaInnerWheel
Transition from outer to inner wheel in EME2.
double m_minBadLArQuality
the minimal cell quality in the LAr for declaring a cell bad
bool m_twoGaussianNoise
if set to true use 2-gaussian noise description for TileCal
ToolHandle< ILArHVFraction > m_larHVFraction

Member Function Documentation

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }

◆ DeclareInterfaceID()

CaloClusterCollectionProcessor::DeclareInterfaceID ( CaloClusterCollectionProcessor ,
1 ,
0  )
inherited

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ execute() [1/2]

StatusCode CaloClusterMomentsMaker::execute ( const EventContext & ctx,
xAOD::CaloClusterContainer * collection ) const
finaloverridevirtual

Execute on an entire collection of clusters.

Parameters
collectionThe container of clusters. param ctx The event context.

Implements CaloClusterCollectionProcessor.

Definition at line 340 of file CaloClusterMomentsMaker.cxx.

343{
344 ATH_MSG_DEBUG("Executing " << name());
345
346 // Maps cell IdentifierHash to cluster index in cluster collection.
347 // Only used when cluster isolation moment is calculated.
348 using clusterIdx_t = std::uint16_t;
349 typedef std::pair<clusterIdx_t, clusterIdx_t> clusterPair_t;
350 std::vector<clusterPair_t> clusterIdx;
351 const clusterIdx_t noCluster = std::numeric_limits<clusterIdx_t>::max();
352
353 const CaloNoise* noise=nullptr;
355 SG::ReadCondHandle<CaloNoise> noiseHdl{m_noiseCDOKey,ctx};
356 noise=*noiseHdl;
357 }
358
359 SG::ReadCondHandle<CaloDetDescrManager> caloMgrHandle{ m_caloMgrKey, ctx };
360 const CaloDetDescrManager* caloDDMgr = *caloMgrHandle;
361 // Counters for number of empty and non-empty neighbor cells per sampling
362 // layer Only used when cluster isolation moment is calculated.
363 int nbEmpty[CaloCell_ID::Unknown];
364 int nbNonEmpty[CaloCell_ID::Unknown];
365
366 // prepare stuff from entire collection in case isolation moment
367 // should be calculated
368
369 if ( m_calculateIsolation ) {
370
371 if (theClusColl->size() >= noCluster) {
372 msg(MSG::ERROR) << "Too many clusters" << endmsg;
373 return StatusCode::FAILURE;
374 }
375
376 // initialize with "empty" values
377 clusterIdx.resize(m_calo_id->calo_cell_hash_max(),
378 clusterPair_t(noCluster, noCluster));
379
380 int iClus = 0;
381 for (xAOD::CaloCluster* theCluster : *theClusColl) {
382 // loop over all cell members and fill cell vector for used cells
383 xAOD::CaloCluster::cell_iterator cellIter = theCluster->cell_begin();
384 xAOD::CaloCluster::cell_iterator cellIterEnd = theCluster->cell_end();
385 for(; cellIter != cellIterEnd; cellIter++ ){
386 CxxUtils::prefetchNext(cellIter, cellIterEnd);
387 const CaloCell* pCell = *cellIter;
388
389 Identifier myId = pCell->ID();
390 IdentifierHash myHashId = m_calo_id->calo_cell_hash(myId);
391 if ( clusterIdx[(unsigned int)myHashId].first != noCluster) {
392 // check weight and assign to current cluster if weight is > 0.5
393 double weight = cellIter.weight();
394 if ( weight > 0.5 )
395 clusterIdx[(unsigned int)myHashId].first = iClus;
396 }
397 else {
398 clusterIdx[(unsigned int)myHashId].first = iClus;
399 }
400 }
401 ++iClus;
402 }
403 }
404
405 // Move allocation of temporary arrays outside the cluster loop.
406 // That way, we don't need to delete and reallocate them
407 // each time through the loop.
408
409 std::vector<CaloClusterMomentsMaker_detail::cellinfo> cellinfo;
410 std::vector<double> maxSampE (CaloCell_ID::Unknown);
411 std::vector<double> myMoments(m_validMoments.size(),0);
412 std::vector<double> myNorms(m_validMoments.size(),0);
413 std::vector<std::tuple<int,int> > nCellsSamp; nCellsSamp.reserve(CaloCell_ID::Unknown);
414 std::vector<IdentifierHash> theNeighbors;
415 // loop over individual clusters
416 xAOD::CaloClusterContainer::iterator clusIter = theClusColl->begin();
417 xAOD::CaloClusterContainer::iterator clusIterEnd = theClusColl->end();
418 int iClus = 0;
419 for( ;clusIter!=clusIterEnd;++clusIter,++iClus) {
420 xAOD::CaloCluster * theCluster = *clusIter;
421
422 double w(0),xc(0),yc(0),zc(0),mx(0),my(0),mz(0),mass(0);
423 double eBad(0),ebad_dac(0),ePos(0),eBadLArQ(0),sumSig2(0),maxAbsSig(0);
424 double eLAr2(0),eLAr2Q(0);
425 double eTile2(0),eTile2Q(0);
426 double eBadLArHV(0);
427 int nbad(0),nbad_dac(0),nBadLArHV(0);
428 unsigned int i,nSigSampl(0);
429 unsigned int theNumOfCells = theCluster->size();
430
431 // these two are needed for the LATERAL moment
432 int iCellMax(-1);
433 int iCellScndMax(-1);
434
435 cellinfo.clear();
436 if (cellinfo.capacity() == 0)
437 cellinfo.reserve (theNumOfCells*2);
438
439 for(i=0;i<(unsigned int)CaloCell_ID::Unknown;i++)
440 maxSampE[i] = 0;
441
442 if ( !m_momentsNames.empty() ) {
443 std::fill (myMoments.begin(), myMoments.end(), 0);
444 std::fill (myNorms.begin(), myNorms.end(), 0);
445 if ( m_calculateIsolation ) {
446 std::fill_n(nbNonEmpty, CaloCell_ID::Unknown, 0);
447 std::fill_n(nbEmpty, CaloCell_ID::Unknown, 0);
448 }
449
450 // loop over all cell members and calculate the center of mass
451 xAOD::CaloCluster::cell_iterator cellIter = theCluster->cell_begin();
452 xAOD::CaloCluster::cell_iterator cellIterEnd = theCluster->cell_end();
453 for(; cellIter != cellIterEnd; cellIter++ ){
454 CaloPrefetch::nextDDE(cellIter, cellIterEnd);
455
456 const CaloCell* pCell = (*cellIter);
457 Identifier myId = pCell->ID();
458 const CaloDetDescrElement* myCDDE = pCell->caloDDE();
459 double ene = pCell->e();
460 if(m_absOpt) ene = std::abs(ene);
461 double weight = cellIter.weight();//theCluster->getCellWeight(cellIter);
462 if ( pCell->badcell() ) {
463 eBad += ene*weight;
464 nbad++;
465 if(ene!=0){
466 ebad_dac+=ene*weight;
467 nbad_dac++;
468 }
469 }
470 else {
471 if ( myCDDE && ! (myCDDE->is_tile())
472 && ((pCell->provenance() & 0x2000) == 0x2000)
473 && !((pCell->provenance() & 0x0800) == 0x0800)) {
474 if ( pCell->quality() > m_minBadLArQuality ) {
475 eBadLArQ += ene*weight;
476 }
477 eLAr2 += ene*weight*ene*weight;
478 eLAr2Q += ene*weight*ene*weight*pCell->quality();
479 }
480 if ( myCDDE && myCDDE->is_tile() ) {
481 uint16_t tq = pCell->quality();
482 uint8_t tq1 = (0xFF00&tq)>>8; // quality in channel 1
483 uint8_t tq2 = (0xFF&tq); // quality in channel 2
484 // reject cells with either 0xFF00 or 0xFF
485 if ( ((tq1&0xFF) != 0xFF) && ((tq2&0xFF) != 0xFF) ) {
486 eTile2 += ene*weight*ene*weight;
487 // take the worse of both qualities (one might be 0 in
488 // 1-channel cases)
489 eTile2Q += ene*weight*ene*weight*(tq1>tq2?tq1:tq2);
490 }
491 }
492 }
493 if ( ene > 0 ) {
494 ePos += ene*weight;
495 }
496
498 const float sigma = m_twoGaussianNoise ?\
499 noise->getEffectiveSigma(pCell->ID(),pCell->gain(),pCell->energy()) : \
500 noise->getNoise(pCell->ID(),pCell->gain());
501
502 sumSig2 += sigma*sigma;
503 // use geometry weighted energy of cell for leading cell significance
504 double Sig = (sigma>0?ene*weight/sigma:0);
505 if (m_useGPUCriteria) {
506 unsigned int thisSampl = myCDDE->getSampling();
507 if ( ( std::abs(Sig) > std::abs(maxAbsSig) ) ||
508 ( std::abs(Sig) == std::abs(maxAbsSig) && thisSampl > nSigSampl ) ||
509 ( std::abs(Sig) == std::abs(maxAbsSig) && thisSampl == nSigSampl && Sig > maxAbsSig ) ) {
510 maxAbsSig = Sig;
511 nSigSampl = thisSampl;
512 }
513
514 }
515 else {
516 if ( std::abs(Sig) > std::abs(maxAbsSig) ) {
517 maxAbsSig = Sig;
518 nSigSampl = myCDDE->getSampling();
519 }
520 }
521 }
522 if ( m_calculateIsolation ) {
523 // get all 2D Neighbours if the cell is not inside another cluster with
524 // larger weight
525
526 IdentifierHash myHashId = m_calo_id->calo_cell_hash(myId);
527 if ( clusterIdx[myHashId].first == iClus ) {
528 theNeighbors.clear();
529 m_calo_id->get_neighbours(myHashId, LArNeighbours::all2D, theNeighbors);
530 for (const auto& nhash: theNeighbors) {
531 clusterPair_t& idx = clusterIdx[nhash];
532
533 // only need to look at each cell once per cluster
534 if ( idx.second == iClus ) continue;
535 idx.second = iClus;
536
537 if ( idx.first == noCluster ) {
538 ++ nbEmpty[m_calo_id->calo_sample(nhash)];
539 } else if ( idx.first != iClus ) {
540 ++ nbNonEmpty[m_calo_id->calo_sample(nhash)];
541 }
542
543 }
544 }
545 }
546
547 if ( myCDDE != nullptr ) {
548 if ( m_nCellsPerSampling ) {
549 CaloCell_ID::CaloSample sam = myCDDE->getSampling();
550 size_t idx((size_t)sam);
551 if ( idx >= nCellsSamp.size() ) { nCellsSamp.resize(idx+1, { 0, 0 } ); }
552 std::get<0>(nCellsSamp[idx])++;
553 // special count for inner wheel cells in EME2
554 if ( sam == CaloCell_ID::EME2 && std::abs(myCDDE->eta()) > m_etaInnerWheel ) { std::get<1>(nCellsSamp[idx])++; }
555 }
556 if ( ene > 0. && weight > 0) {
557 // get all geometric information needed ...
558 cellinfo.push_back (CaloClusterMomentsMaker_detail::cellinfo {
559 .x = myCDDE->x(),
560 .y = myCDDE->y(),
561 .z = myCDDE->z(),
562 .energy = ene*weight,
563 .eta = myCDDE->eta(),
564 .phi = myCDDE->phi(),
565 .r = 0, // These two are filled in later.
566 .lambda = 0,
567 .volume = myCDDE->volume(),
568 .sample = myCDDE->getSampling(),
569 .identifier = cellIter.index()
570 //Using the index instead of the hash ID as disambiguation criterion
571 //is relevant for the GPU now that we no longer assume the two are the same,
572 //and this gives us better performance.
573 });
574
575 CaloClusterMomentsMaker_detail::cellinfo& ci = cellinfo.back();
576
577 if ( ci.energy > maxSampE[(unsigned int)ci.sample] )
578 maxSampE[(unsigned int)ci.sample] = ci.energy;
579
580 if (m_useGPUCriteria) {
581 if (iCellMax < 0 ||
582 ci.energy > cellinfo[iCellMax].energy ||
583 (ci.energy == cellinfo[iCellMax].energy && ci.identifier > cellinfo[iCellMax].identifier) ) {
584 iCellScndMax = iCellMax;
585 iCellMax = cellinfo.size()-1;
586 }
587 else if (iCellScndMax < 0 ||
588 ci.energy > cellinfo[iCellScndMax].energy ||
589 (ci.energy == cellinfo[iCellScndMax].energy && ci.identifier > cellinfo[iCellScndMax].identifier) )
590 {
591 iCellScndMax = cellinfo.size()-1;
592 }
593 }
594 else {
595 if (iCellMax < 0 || ci.energy > cellinfo[iCellMax].energy ) {
596 iCellScndMax = iCellMax;
597 iCellMax = cellinfo.size()-1;
598 }
599 else if (iCellScndMax < 0 ||
600 ci.energy > cellinfo[iCellScndMax].energy )
601 {
602 iCellScndMax = cellinfo.size()-1;
603 }
604 }
605
606 xc += ci.energy*ci.x;
607 yc += ci.energy*ci.y;
608 zc += ci.energy*ci.z;
609
610 double dir = ci.x*ci.x+ci.y*ci.y+ci.z*ci.z;
611
612 if ( dir > 0) {
613 dir = sqrt(dir);
614 dir = 1./dir;
615 }
616 mx += ci.energy*ci.x*dir;
617 my += ci.energy*ci.y*dir;
618 mz += ci.energy*ci.z*dir;
619
620 w += ci.energy;
621 } // cell has E>0 and weight != 0
622 } // cell has valid DDE
623 } //end of loop over all cells
625 const auto hvFrac=m_larHVFraction->getLArHVFrac(theCluster->getCellLinks(),ctx);
626 eBadLArHV= hvFrac.first;
627 nBadLArHV=hvFrac.second;
628 }
629
630 if ( w > 0 || (m_useGPUCriteria && w >= 0) ) {
631 mass = w*w - mx*mx - my*my - mz*mz;
632 if ( mass > 0) {
633 mass = sqrt(mass);
634 }
635 else {
636 // make mass negative if m^2 was negative
637 mass = -sqrt(-mass);
638 }
639
640 if (w == 0) {
641 w = 1.0;
642 }
643
644 xc/=w;
645 yc/=w;
646 zc/=w;
647 Amg::Vector3D showerCenter(xc,yc,zc);
648 w=0;
649
650
651 //log << MSG::WARNING << "Found bad cells " << xbad_dac << " " << ybad_dac << " " << zbad_dac << " " << ebad_dac << endmsg;
652 //log << MSG::WARNING << "Found Cluster " << xbad_dac << " " << ybad_dac << " " << zbad_dac << " " << endmsg;
653 // shower axis is just the vector pointing from the IP to the shower center
654 // in case there are less than 3 cells in the cluster
655
656 Amg::Vector3D showerAxis(xc,yc,zc);
657 Amg::setMag(showerAxis,1.0);
658
659 // otherwise the principal direction with the largest absolute
660 // eigenvalue will be used unless it's angle w.r.t. the vector pointing
661 // from the IP to the shower center is larger than allowed by the
662 // property m_maxAxisAngle
663
664 double angle(0),deltaPhi(0),deltaTheta(0);
665 if ( cellinfo.size() > 2 ) {
666 Eigen::Matrix3d C=Eigen::Matrix3d::Zero();
667 for(const CaloClusterMomentsMaker_detail::cellinfo& ci : cellinfo) {
668 const double e2 = ci.energy * ci.energy;
669
670 C(0,0) += e2*(ci.x-xc)*(ci.x-xc);
671 C(1,0) += e2*(ci.x-xc)*(ci.y-yc);
672 C(2,0) += e2*(ci.x-xc)*(ci.z-zc);
673
674 C(1,1) += e2*(ci.y-yc)*(ci.y-yc);
675 C(2,1) += e2*(ci.y-yc)*(ci.z-zc);
676
677 C(2,2) += e2*(ci.z-zc)*(ci.z-zc);
678 w += e2;
679 }
680 C/=(w != 0 ? w : 1.0);
681
682 Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> eigensolver(C);
683 if (eigensolver.info() != Eigen::Success) {
684 msg(MSG::WARNING) << "Failed to compute Eigenvalues -> Can't determine shower axis" << endmsg;
685 }
686 else {
687 // don't use the principal axes if at least one of the 3
688 // diagonal elements is 0
689
690 const Eigen::Vector3d& S=eigensolver.eigenvalues();
691 const Eigen::Matrix3d& U=eigensolver.eigenvectors();
692
693 const double epsilon = 1.E-6;
694
695 if ( std::abs(S[0]) >= epsilon && std::abs(S[1]) >= epsilon && std::abs(S[2]) >= epsilon ) {
696
697 Amg::Vector3D prAxis(showerAxis);
698 int iEigen = -1;
699
700 for (i=0;i<3;i++) {
701 Amg::Vector3D tmpAxis=U.col(i);
702
703 // calculate the angle
704 double tmpAngle=Amg::angle(tmpAxis,showerAxis);
705
706 if ( tmpAngle > 90*deg ) {
707 tmpAngle = 180*deg - tmpAngle;
708 tmpAxis = -tmpAxis;
709 }
710
711 if ( iEigen == -1 || tmpAngle < angle ) {
712 iEigen = i;
713 angle = tmpAngle;
714 prAxis = tmpAxis;
715 }
716 }//end for loop
717
718 // calculate theta and phi angle differences
719
720 deltaPhi = CaloPhiRange::diff(showerAxis.phi(),prAxis.phi());
721
722 deltaTheta = showerAxis.theta() - prAxis.theta();
723
724 // check the angle
725
726 if ( angle < m_maxAxisAngle ) {
727 showerAxis = prAxis;
728 }
729 else
730 ATH_MSG_DEBUG("principal Direction (" << prAxis[Amg::x] << ", "
731 << prAxis[Amg::y] << ", " << prAxis[Amg::z] << ") deviates more than "
732 << m_maxAxisAngle*(1./deg)
733 << " deg from IP-to-ClusterCenter-axis (" << showerAxis[Amg::x] << ", "
734 << showerAxis[Amg::y] << ", " << showerAxis[Amg::z] << ")");
735 }//end if std::abs(S)<epsilon
736 else {
737 ATH_MSG_DEBUG("Eigenvalues close to 0, do not use principal axis");
738 }
739 }//end got eigenvalues
740 } //end if cellinfo.size()>2
741
742 ATH_MSG_DEBUG("Shower Axis = (" << showerAxis[Amg::x] << ", "
743 << showerAxis[Amg::y] << ", " << showerAxis[Amg::z] << ")");
744
745
746 // calculate radial distance from and the longitudinal distance
747 // along the shower axis for each cell. The cluster center is
748 // at r=0 and lambda=0
749
750 for (auto& ci : cellinfo) {
751 const Amg::Vector3D currentCell(ci.x,ci.y,ci.z);
752 // calculate distance from shower axis r
753 ci.r = ((currentCell-showerCenter).cross(showerAxis)).mag();
754 // calculate distance from shower center along shower axis
755 ci.lambda = (currentCell-showerCenter).dot(showerAxis);
756 }
757
758 // loop over all positive energy cells and calculate all desired moments
759
760 // define common norm for all simple moments
761 double commonNorm = 0;
762 double phi0 = cellinfo.size() > 0 ? cellinfo[0].phi : 0;
763
764 for(unsigned i=0;i<cellinfo.size();i++) {
765 const CaloClusterMomentsMaker_detail::cellinfo& ci = cellinfo[i];
766 // loop over all valid moments
767 commonNorm += ci.energy;
768 for(size_t iMoment = 0, size = m_validMoments.size();
769 iMoment != size;
770 ++ iMoment)
771 {
772 // now calculate the actual moments
773 switch (m_validMoments[iMoment]) {
774 case FIRST_ETA:
775 myMoments[iMoment] += ci.energy*ci.eta;
776 break;
777 case FIRST_PHI:
778 // first cell decides the sign in order to avoid
779 // overlap problem at phi = -pi == +pi
780 // need to be normalized to the range [-pi,+pi] in the end
781 myMoments[iMoment] += ci.energy * proxim (ci.phi, phi0);
782 break;
783 case SECOND_R:
784 myMoments[iMoment] += ci.energy*ci.r*ci.r;
785 break;
786 case SECOND_LAMBDA:
787 myMoments[iMoment] += ci.energy*ci.lambda*ci.lambda;
788 break;
789 case LATERAL:
790 if ( (int)i != iCellMax && (int)i != iCellScndMax ) {
791 myMoments[iMoment] += ci.energy*ci.r*ci.r;
792 myNorms[iMoment] += ci.energy*ci.r*ci.r;
793 }
794 else {
795 double rm = ci.r;
796 if ( rm < m_minRLateral )
797 rm = m_minRLateral;
798 myNorms[iMoment] += rm*rm*ci.energy;
799 }
800 break;
801 case LONGITUDINAL:
802 if ( (int)i != iCellMax && (int)i != iCellScndMax ) {
803 myMoments[iMoment] += ci.energy*ci.lambda*ci.lambda;
804 myNorms[iMoment] += ci.energy*ci.lambda*ci.lambda;
805 }
806 else {
807 double lm = ci.lambda;
808 if ( lm < m_minLLongitudinal )
810 myNorms[iMoment] += lm*lm*ci.energy;
811 }
812 break;
813 case FIRST_ENG_DENS:
814 if ( ci.volume > 0 ) {
815 myMoments[iMoment] += ci.energy*ci.energy/ci.volume;
816 myNorms[iMoment] += ci.energy;
817 }
818 break;
819 case SECOND_ENG_DENS:
820 if ( ci.volume > 0 ) {
821 myMoments[iMoment] += ci.energy*std::pow(ci.energy/ci.volume,2);
822 myNorms[iMoment] += ci.energy;
823 }
824 break;
825 case ENG_FRAC_EM:
826 if ( ci.sample == CaloCell_ID::EMB1
827 || ci.sample == CaloCell_ID::EMB2
828 || ci.sample == CaloCell_ID::EMB3
829 || ci.sample == CaloCell_ID::EME1
830 || ci.sample == CaloCell_ID::EME2
831 || ci.sample == CaloCell_ID::EME3
832 || ci.sample == CaloCell_ID::FCAL0 )
833 myMoments[iMoment] += ci.energy;
834 break;
835 case ENG_FRAC_MAX:
836 if ( (int)i == iCellMax )
837 myMoments[iMoment] = ci.energy;
838 break;
839 case PTD:
840 // do not convert to pT since clusters are small and
841 // there is virtually no difference and cosh just costs
842 // time ...
843 myMoments[iMoment] += ci.energy*ci.energy;
844 myNorms[iMoment] += ci.energy;
845 break;
846 default:
847 // nothing to be done for other moments
848 break;
849 }
850 }
851 } //end of loop over cell
852
853
854 // assign moments which don't need the loop over the cells
855 for (size_t iMoment = 0, size = m_validMoments.size();
856 iMoment != size;
857 ++ iMoment)
858 {
859 // now calculate the actual moments
860 switch (m_validMoments[iMoment]) {
861 case FIRST_ETA:
862 case FIRST_PHI:
863 case SECOND_R:
864 case SECOND_LAMBDA:
865 case ENG_FRAC_EM:
866 case ENG_FRAC_MAX:
867 myNorms[iMoment] = commonNorm;
868 break;
869 case DELTA_PHI:
870 myMoments[iMoment] = deltaPhi;
871 break;
872 case DELTA_THETA:
873 myMoments[iMoment] = deltaTheta;
874 break;
875 case DELTA_ALPHA:
876 myMoments[iMoment] = angle;
877 break;
878 case CENTER_X:
879 myMoments[iMoment] = showerCenter.x();
880 break;
881 case CENTER_Y:
882 myMoments[iMoment] = showerCenter.y();
883 break;
884 case CENTER_Z:
885 myMoments[iMoment] = showerCenter.z();
886 break;
887 case CENTER_MAG:
888 myMoments[iMoment] = showerCenter.mag();
889 break;
890 case CENTER_LAMBDA:
891 // calculate the longitudinal distance along the shower axis
892 // of the shower center from the calorimeter start
893
894 // first need calo boundary at given eta phi try LAREM barrel
895 // first, then LAREM endcap OW, then LAREM endcap IW, then
896 // FCal
897 {
898 double r_calo(0),z_calo(0),lambda_c(0);
899 r_calo = m_caloDepthTool->get_entrance_radius(CaloCell_ID::EMB1,
900 showerCenter.eta(),
901 showerCenter.phi(),
902 caloDDMgr);
903 if ( r_calo == 0 ) {
904 z_calo = m_caloDepthTool->get_entrance_z(CaloCell_ID::EME1,
905 showerCenter.eta(),
906 showerCenter.phi(),
907 caloDDMgr);
908 if ( z_calo == 0 )
909 z_calo = m_caloDepthTool->get_entrance_z(CaloCell_ID::EME2,
910 showerCenter.eta(),
911 showerCenter.phi(),
912 caloDDMgr);
913 if ( z_calo == 0 )
914 z_calo = m_caloDepthTool->get_entrance_z(CaloCell_ID::FCAL0,
915 showerCenter.eta(),
916 showerCenter.phi(),
917 caloDDMgr);
918 if ( z_calo == 0 ) // for H6 TB without EMEC outer wheel
919 z_calo = m_caloDepthTool->get_entrance_z(CaloCell_ID::HEC0,
920 showerCenter.eta(),
921 showerCenter.phi(),
922 caloDDMgr);
923 if ( z_calo != 0 && showerAxis.z() != 0 ) {
924 lambda_c = std::abs((z_calo-showerCenter.z())/showerAxis.z());
925 }
926 }
927 else {
928 double r_s2 = showerAxis.x()*showerAxis.x()
929 +showerAxis.y()*showerAxis.y();
930 double r_cs = showerAxis.x()*showerCenter.x()
931 +showerAxis.y()*showerCenter.y();
932 double r_cr = showerCenter.x()*showerCenter.x()
933 +showerCenter.y()*showerCenter.y()-r_calo*r_calo;
934 if ( r_s2 > 0 ) {
935 double det = r_cs*r_cs/(r_s2*r_s2) - r_cr/r_s2;
936 if ( det > 0 ) {
937 det = sqrt(det);
938 double l1(-r_cs/r_s2);
939 double l2(l1);
940 l1 += det;
941 l2 -= det;
942 if ( std::abs(l1) < std::abs(l2) )
943 lambda_c = std::abs(l1);
944 else
945 lambda_c = std::abs(l2);
946 }
947 }
948 }
949 myMoments[iMoment] = lambda_c;
950 }
951 break;
952 case ENG_FRAC_CORE:
953 for(i=0;i<(int)CaloCell_ID::Unknown;i++)
954 myMoments[iMoment] += maxSampE[i];
955 myNorms[iMoment] = commonNorm;
956 break;
957 case ISOLATION:
958 {
959 // loop over empty and filled perimeter cells and
960 // get a weighted ratio by means of energy fraction per layer
961 for(unsigned int i=0; i != CaloSampling::Unknown; ++ i) {
963 if (theCluster->hasSampling(s)) {
964 const double eSample = theCluster->eSample(s);
965 if (eSample > 0) {
966 int nAll = nbEmpty[i]+nbNonEmpty[i];
967 if (nAll > 0) {
968 myMoments[iMoment] += (eSample*nbEmpty[i])/nAll;
969 myNorms[iMoment] += eSample;
970 }
971 }//end of eSample>0
972 }//end has sampling
973 }//end loop over samplings
974 }
975 break;
976 case ENG_BAD_CELLS:
977 myMoments[iMoment] = eBad;
978 break;
979 case N_BAD_CELLS:
980 myMoments[iMoment] = nbad;
981 break;
982 case N_BAD_CELLS_CORR:
983 myMoments[iMoment] = nbad_dac;
984 break;
985 case BAD_CELLS_CORR_E:
986 myMoments[iMoment] = ebad_dac;
987 break;
988 case BADLARQ_FRAC:
989 myMoments[iMoment] = eBadLArQ/(theCluster->e()!=0.?theCluster->e():1.);
990 break;
991 case ENG_POS:
992 myMoments[iMoment] = ePos;
993 break;
994 case SIGNIFICANCE:
995 myMoments[iMoment] = (sumSig2>0?theCluster->e()/sqrt(sumSig2):0.);
996 break;
997 case CELL_SIGNIFICANCE:
998 myMoments[iMoment] = maxAbsSig;
999 break;
1000 case CELL_SIG_SAMPLING:
1001 myMoments[iMoment] = nSigSampl;
1002 break;
1003 case AVG_LAR_Q:
1004 myMoments[iMoment] = eLAr2Q/(eLAr2>0?eLAr2:1);
1005 break;
1006 case AVG_TILE_Q:
1007 myMoments[iMoment] = eTile2Q/(eTile2>0?eTile2:1);
1008 break;
1009 case ENG_BAD_HV_CELLS:
1010 myMoments[iMoment] = eBadLArHV;
1011 break;
1012 case N_BAD_HV_CELLS:
1013 myMoments[iMoment] = nBadLArHV;
1014 break;
1015 case PTD:
1016 myMoments[iMoment] = sqrt(myMoments[iMoment]);
1017 break;
1018 case MASS:
1019 myMoments[iMoment] = mass;
1020 break;
1021 default:
1022 // nothing to be done for other moments
1023 break;
1024 }
1025 }
1026 }
1027
1028 // normalize moments and copy to Cluster Moment Store
1029 size_t size= m_validMoments.size();
1030 for (size_t iMoment = 0; iMoment != size; ++iMoment) {
1032 if ( myNorms[iMoment] != 0 )
1033 myMoments[iMoment] /= myNorms[iMoment];
1034 if ( moment == FIRST_PHI )
1035 myMoments[iMoment] = CaloPhiRange::fix(myMoments[iMoment]);
1036 theCluster->insertMoment(moment,myMoments[iMoment]);
1037 } // loop on moments for cluster
1038 } // check on requested moments
1039 // check on second moment of time if requested
1040 if ( m_secondTime ) { theCluster->insertMoment(SECOND_TIME,theCluster->secondTime()); }
1041 // check on number of cells per sampling moment if requested
1042 if ( m_nCellsPerSampling ) {
1043 for ( size_t isam(0); isam < nCellsSamp.size(); ++isam ) {
1044 theCluster->setNumberCellsInSampling((CaloCell_ID::CaloSample)isam,std::get<0>(nCellsSamp.at(isam)),false);
1045 if ( isam == (size_t)CaloCell_ID::EME2 && std::get<1>(nCellsSamp.at(isam)) > 0 ) {
1046 theCluster->setNumberCellsInSampling((CaloCell_ID::CaloSample)isam,std::get<1>(nCellsSamp.at(isam)),true);
1047 }
1048 } // loop on samplings
1049 nCellsSamp.clear();
1050 }
1051 } // loop on clusters
1052
1053 return StatusCode::SUCCESS;
1054}
Scalar deltaPhi(const MatrixBase< Derived > &vec) const
#define endmsg
#define ATH_MSG_DEBUG(x)
double angle(const GeoTrf::Vector2D &a, const GeoTrf::Vector2D &b)
MsgStream & msg() const
CaloSampling::CaloSample CaloSample
Definition CaloCell_ID.h:53
virtual double e() const override final
get energy (data member) (synonym to method energy()
Definition CaloCell.h:333
double energy() const
get energy (data member)
Definition CaloCell.h:327
const CaloDetDescrElement * caloDDE() const
get pointer to CaloDetDescrElement (data member)
Definition CaloCell.h:321
uint16_t provenance() const
get provenance (data member)
Definition CaloCell.h:354
uint16_t quality() const
get quality (data member)
Definition CaloCell.h:348
virtual bool badcell() const
check is cell is dead
Definition CaloCell.cxx:144
CaloGain::CaloGain gain() const
get gain (data member )
Definition CaloCell.h:361
Identifier ID() const
get ID (from cached data member) non-virtual and inline for fast access
Definition CaloCell.h:295
SG::ReadCondHandleKey< CaloDetDescrManager > m_caloMgrKey
bool m_secondTime
Retrieve second moment of cell times and store as moment.
std::vector< xAOD::CaloCluster::MomentType > m_validMoments
set of moments which will be calculated.
Gaudi::Property< bool > m_useGPUCriteria
bool m_nCellsPerSampling
store number of cells per sampling layer as moment
SG::ReadCondHandleKey< CaloNoise > m_noiseCDOKey
Key of the CaloNoise Conditions data object.
CaloCell_ID::CaloSample getSampling() const
cell sampling
static double fix(double phi)
static double diff(double phi1, double phi2)
simple phi1 - phi2 calculation, but result is fixed to respect range.
DataModel_detail::iterator< DataVector > iterator
Definition DataVector.h:842
const CaloClusterCellLink * getCellLinks() const
Get a pointer to the CaloClusterCellLink object (const version)
size_t size() const
size method (forwarded from CaloClusterCellLink obj)
CaloClusterCellLink::iterator cell_iterator
Iterator of the underlying CaloClusterCellLink (non-const version)
virtual double e() const
The total energy of the particle.
void insertMoment(MomentType type, double value)
const_cell_iterator cell_end() const
float eSample(const CaloSample sampling) const
flt_t secondTime() const
Access second moment of cell timing distribution.
MomentType
Enums to identify different moments.
const_cell_iterator cell_begin() const
Iterator of the underlying CaloClusterCellLink (const version)
void setNumberCellsInSampling(CaloSampling::CaloSample samp, int ncells, bool isInnerWheel=false)
Set the number of cells in a sampling layer.
CaloSampling::CaloSample CaloSample
bool hasSampling(const CaloSample s) const
Checks if certain smapling contributes to cluster.
struct color C
void setMag(Amg::Vector3D &v, double mag)
scales the vector length without changing the angles
double angle(const Amg::Vector3D &v1, const Amg::Vector3D &v2)
calculates the opening angle between two vectors
Eigen::Matrix< double, 3, 1 > Vector3D
void nextDDE(Iter iter, Iter endIter)
Prefetch next CaloDDE.
void prefetchNext(Iter iter, Iter endIter)
Prefetch next object in sequence.
Definition prefetch.h:130
bool FIRST_ENG_DENS(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool CENTER_MAG(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool CENTER_LAMBDA(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool SECOND_R(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool SECOND_LAMBDA(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
dot(G, fn, nodesToHighlight=[])
Definition dot.py:5
double e2(const xAOD::CaloCluster &cluster)
return the uncorrected cluster energy in 2nd sampling
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
setWord1 uint16_t
double proxim(double b, double a)
Definition proxim.h:17

◆ execute() [2/2]

virtual StatusCode CaloClusterCollectionProcessor::execute ( xAOD::CaloClusterContainer * collection)
inlinefinalvirtual

Execute on an entire collection of clusters.

Parameters
collectionThe container of clusters. (deprecated)

Reimplemented from CaloClusterCollectionProcessor.

Definition at line 50 of file CaloClusterCollectionProcessor.h.

51 {
52 return execute (Gaudi::Hive::currentContext(), collection);
53 }
virtual StatusCode execute(const EventContext &ctx, xAOD::CaloClusterContainer *theClusColl) const override final
Execute on an entire collection of clusters.

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ finalize()

StatusCode CaloClusterMomentsMaker::finalize ( )
overridevirtual

Definition at line 313 of file CaloClusterMomentsMaker.cxx.

314{
315 return StatusCode::SUCCESS;
316}

◆ initialize()

StatusCode CaloClusterMomentsMaker::initialize ( )
overridevirtual

Definition at line 183 of file CaloClusterMomentsMaker.cxx.

184{
185 xAOD::CaloCluster dummyCluster;
186
187 // loop list of requested moments
188 std::string::size_type nstr(0); int nmom(0);
189 for (const auto& mom : m_momentsNames) {
190 ATH_MSG_DEBUG("Moment " << mom << " requested");
191 // check if moment is known (enumerator available)
192 auto fmap(momentNameToEnumMap.find(mom));
193 if (fmap != momentNameToEnumMap.end()) {
194 // valid moment found
195 nstr = std::max(nstr, mom.length());
196 ++nmom;
197 if (fmap->second == SECOND_TIME) {
198 // special flag for second moment of cell times - this moment is not
199 // calculated in this tool! Do not add to internal (!) valid moments
200 // list. Its value is available from xAOD::CaloCluster::secondTime()!
201 m_secondTime = true;
202
203 // Make sure the variable used for the moment is declared
204 // to the auxiliary variable registry. Otherwise, if we don't
205 // set the moment for the first event (perhaps because there
206 // are no clusters), then we can get warnings from AuxSelection.
207 (void)dummyCluster.getMomentValue (fmap->second);
208 } else if (fmap->second == NCELL_SAMPLING) {
209 // flag indicates if number of cells in a sampling should be counted.
210 // This is a vector of integers counts that is filled in this tool but
211 // does not need any post-processing (e.g. normalization). It is not
212 // added to the valid moments list for this reason.
213 ATH_MSG_DEBUG("moment " << fmap->first << " found");
214 m_nCellsPerSampling = true;
216
217 // Make sure the variable used for the moment is declared
218 // to the auxiliary variable registry.
219 (void)dummyCluster.retrieveMoment (fmap->second, cellsdum);
220 } else if (fmap->second == EM_PROBABILITY) {
222 << " not calculated in this tool - misconfiguration?");
223 } else {
224 // Make sure the variable used for the moment is declared
225 // to the auxiliary variable registry.
226 (void)dummyCluster.getMomentValue (fmap->second);
227
228 // all other valid moments
229 m_validMoments.push_back(fmap->second);
230 // flag some special requests
231 switch (fmap->second) {
232 case SIGNIFICANCE:
233 case CELL_SIGNIFICANCE:
235 break;
236 case ISOLATION:
238 break;
239 case ENG_BAD_HV_CELLS:
241 break;
242 default:
243 break;
244 } // set special processing flags
245 } // moment calculated with this tool
246 } else {
247 ATH_MSG_ERROR("Moment name " << mom << " not known; known moments are:");
248 char buffer[128];
249 std::string::size_type lstr(nstr);
250 // determine field size
251 for (const auto& fmom : momentNameToEnumMap) {
252 lstr = std::max(lstr, fmom.first.length());
253 }
254 // print available moments
255 for (const auto& fmom : momentNameToEnumMap) {
256 sprintf(buffer, "moment name: %-*.*s - enumerator: %i", (int)lstr,
257 (int)lstr, fmom.first.c_str(), (int)fmom.second);
258 ATH_MSG_ERROR(buffer);
259 }
260 auto fmom(momentNameToEnumMap.find("SECOND_TIME"));
261 sprintf(buffer, "moment name: %-*.*s - enumerator: %i", (int)nstr,
262 (int)nstr, fmom->first.c_str(), (int)fmom->second);
263 ATH_MSG_ERROR(buffer);
264 fmom = momentNameToEnumMap.find("NCELL_SAMPLING");
265 sprintf(buffer, "moment name: %-*.*s - enumerator: %i", (int)nstr,
266 (int)nstr, fmom->first.c_str(), (int)fmom->second);
267 ATH_MSG_ERROR(buffer);
268 return StatusCode::FAILURE;
269 } // found unknown moment name
270 } // loop configured moment names
271
272 // sort and remove duplicates
274 m_validMoments.erase(
276 m_validMoments.end());
277
278 // print configured moments
279 ATH_MSG_INFO("Construct and save " << nmom << " cluster moments: ");
280 char buffer[128];
281 for (auto menum : m_validMoments) {
282 sprintf(buffer, "moment name: %-*.*s - enumerator: %i", (int)nstr,
283 (int)nstr, momentEnumToNameMap.at(menum).c_str(), (int)menum);
284 ATH_MSG_INFO(buffer);
285 }
286 if (m_secondTime) {
287 auto fmom(momentNameToEnumMap.find("SECOND_TIME"));
288 sprintf(buffer, "moment name: %-*.*s - enumerator: %i (save only)",
289 (int)nstr, (int)nstr, fmom->first.c_str(), (int)fmom->second);
290 ATH_MSG_INFO(buffer);
291 }
293 auto fmom(momentNameToEnumMap.find("NCELL_SAMPLING"));
294 sprintf(buffer, "moment name: %-*.*s - enumerator: %i", (int)nstr,
295 (int)nstr, fmom->first.c_str(), (int)fmom->second);
296 ATH_MSG_INFO(buffer);
297 }
298
299 // retrieve CaloCell ID server
300 CHECK(detStore()->retrieve(m_calo_id,"CaloCell_ID"));
301
302 // retrieve the calo depth tool
303 CHECK(m_caloDepthTool.retrieve());
304 ATH_CHECK(m_caloMgrKey.initialize());
305
306 // retrieve specific servers and tools for selected processes
308 if (m_calculateLArHVFraction) { ATH_CHECK(m_larHVFraction.retrieve()); } else { m_larHVFraction.disable(); }
309
310 return StatusCode::SUCCESS;
311}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define CHECK(...)
Evaluate an expression and check for errors.
const ServiceHandle< StoreGateSvc > & detStore() const
bool retrieveMoment(MomentType type, double &value) const
Retrieve individual moment.
double getMomentValue(MomentType type) const
Retrieve individual moment - no check for existance! Returns -999 on error.
std::vector< uint16_t > ncells_store_t
Store type for number-of-cells-in-sampling counter.
retrieve(aClass, aKey=None)
Definition PyKernel.py:110
DataModel_detail::iterator< DVL > unique(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of unique for DataVector/List.
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ msg()

MsgStream & AthCommonMsg< AlgTool >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msgLvl()

bool AthCommonMsg< AlgTool >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< AlgTool > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sysInitialize()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysInitialize ( )
overridevirtualinherited

Perform system initialization for an algorithm.

We override this to declare all the elements of handle key arrays at the end of initialization. See comments on updateVHKA.

Reimplemented in asg::AsgMetadataTool, AthCheckedComponent< AthAlgTool >, AthCheckedComponent<::AthAlgTool >, and DerivationFramework::CfAthAlgTool.

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

Member Data Documentation

◆ m_absOpt

bool CaloClusterMomentsMaker::m_absOpt
private

if set to true use abs E value of cells to calculate cluster moments

Definition at line 146 of file CaloClusterMomentsMaker.h.

◆ m_calculateIsolation

bool CaloClusterMomentsMaker::m_calculateIsolation
private

Set to true if cluster isolation is to be calculated.

Definition at line 115 of file CaloClusterMomentsMaker.h.

◆ m_calculateLArHVFraction

bool CaloClusterMomentsMaker::m_calculateLArHVFraction
private

Set to true to calculate E and N of cells affected by LAr HV corrections.

Definition at line 118 of file CaloClusterMomentsMaker.h.

◆ m_calculateSignificance

bool CaloClusterMomentsMaker::m_calculateSignificance
private

Set to true if significance moments are need.

Definition at line 112 of file CaloClusterMomentsMaker.h.

◆ m_calo_id

const CaloCell_ID* CaloClusterMomentsMaker::m_calo_id
private

Definition at line 79 of file CaloClusterMomentsMaker.h.

◆ m_caloDepthTool

ToolHandle<CaloDepthTool> CaloClusterMomentsMaker::m_caloDepthTool
private

Definition at line 125 of file CaloClusterMomentsMaker.h.

◆ m_caloMgrKey

SG::ReadCondHandleKey<CaloDetDescrManager> CaloClusterMomentsMaker::m_caloMgrKey
private
Initial value:
{
this,
"CaloDetDescrManager",
"CaloDetDescrManager"
}

Definition at line 127 of file CaloClusterMomentsMaker.h.

127 {
128 this,
129 "CaloDetDescrManager",
130 "CaloDetDescrManager"
131 };

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_etaInnerWheel

double CaloClusterMomentsMaker::m_etaInnerWheel = { 2.52 }
private

Transition from outer to inner wheel in EME2.

Definition at line 158 of file CaloClusterMomentsMaker.h.

158{ 2.52 };

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_larHVFraction

ToolHandle<ILArHVFraction> CaloClusterMomentsMaker::m_larHVFraction
private

Definition at line 138 of file CaloClusterMomentsMaker.h.

◆ m_maxAxisAngle

double CaloClusterMomentsMaker::m_maxAxisAngle
private

the maximal allowed deviation from the IP-to-ClusterCenter-axis.

Definition at line 83 of file CaloClusterMomentsMaker.h.

◆ m_minBadLArQuality

double CaloClusterMomentsMaker::m_minBadLArQuality
private

the minimal cell quality in the LAr for declaring a cell bad

This defines the minimal quality (large values mean worse shape) a cell needs to exceed in order to be considered as not compatible with a normal ionization signal.

Definition at line 109 of file CaloClusterMomentsMaker.h.

◆ m_minLLongitudinal

double CaloClusterMomentsMaker::m_minLLongitudinal
private

the minimal \(\lambda\) in the definition of the Longitudinal moment

This defines the minimal distance along the shower axis from the cluster center the two leading cells might have before this value is used instead of their real distance in the normalization of the LONGITUDINAL moment.

Definition at line 101 of file CaloClusterMomentsMaker.h.

◆ m_minRLateral

double CaloClusterMomentsMaker::m_minRLateral
private

the minimal \(r\) in the definition of the Lateral moment

This defines the minimal distance the two leading cells might have before this value is used instead of their real distance in the normalization of the LATERAL moment.

Definition at line 91 of file CaloClusterMomentsMaker.h.

◆ m_momentsNames

std::vector<std::string> CaloClusterMomentsMaker::m_momentsNames
private

vector holding the input list of names of moments to calculate.

This is the list of desired names of moments given in the jobOptions.

Definition at line 69 of file CaloClusterMomentsMaker.h.

◆ m_momentsNamesAOD

std::string CaloClusterMomentsMaker::m_momentsNamesAOD
private

Not used anymore (with xAOD), but required when configured from COOL.

Definition at line 141 of file CaloClusterMomentsMaker.h.

◆ m_nCellsPerSampling

bool CaloClusterMomentsMaker::m_nCellsPerSampling = { false }
private

store number of cells per sampling layer as moment

Definition at line 154 of file CaloClusterMomentsMaker.h.

154{ false };

◆ m_noiseCDOKey

SG::ReadCondHandleKey<CaloNoise> CaloClusterMomentsMaker::m_noiseCDOKey {this,"CaloNoiseKey","totalNoise","SG Key of CaloNoise data object"}
private

Key of the CaloNoise Conditions data object.

Typical values are '"electronicNoise', 'pileupNoise', or '"totalNoise' (default)

Definition at line 136 of file CaloClusterMomentsMaker.h.

136{this,"CaloNoiseKey","totalNoise","SG Key of CaloNoise data object"};

◆ m_secondTime

bool CaloClusterMomentsMaker::m_secondTime = { false }
private

Retrieve second moment of cell times and store as moment.

Definition at line 150 of file CaloClusterMomentsMaker.h.

150{ false };

◆ m_twoGaussianNoise

bool CaloClusterMomentsMaker::m_twoGaussianNoise
private

if set to true use 2-gaussian noise description for TileCal

Definition at line 123 of file CaloClusterMomentsMaker.h.

◆ m_useGPUCriteria

Gaudi::Property<bool> CaloClusterMomentsMaker::m_useGPUCriteria {this, "UseGPUCriteria", false, "Adopt a set of criteria that is consistent with the GPU implementation."}
private

Definition at line 161 of file CaloClusterMomentsMaker.h.

161{this, "UseGPUCriteria", false, "Adopt a set of criteria that is consistent with the GPU implementation."};

◆ m_validMoments

std::vector<xAOD::CaloCluster::MomentType> CaloClusterMomentsMaker::m_validMoments
private

set of moments which will be calculated.

This set will hold each valid enum indexed if its name was found on the input list (m_momentsNames) and in the list of valid moment names (m_validNames).

Definition at line 77 of file CaloClusterMomentsMaker.h.

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< AlgTool > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< AlgTool > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.


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