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

Concrete tool that calculates calibration hit truth energies in xAOD::CaloCluster. More...

#include <CaloCalibClusterTruthAttributerTool.h>

Inheritance diagram for CaloCalibClusterTruthAttributerTool:
Collaboration diagram for CaloCalibClusterTruthAttributerTool:

Public Member Functions

 CaloCalibClusterTruthAttributerTool (const std::string &type, const std::string &name, const IInterface *parent)
virtual ~CaloCalibClusterTruthAttributerTool ()
virtual StatusCode calculateTruthEnergies (const xAOD::CaloCluster &theCaloCluster, unsigned int numTruthParticles, const std::map< Identifier, std::vector< const CaloCalibrationHit * > > &identifierToCaloHitMap, std::vector< std::pair< unsigned int, double > > &barcodeTrueCalHitEnergy) const override
 This calculates the truth energies of N leading truth particles in a topocluster.

Private Attributes

Gaudi::Property< bool > m_fullTruthEnergy {this,"storeFullTruthEnergy",false,"Toggle storage of invisible and escaped energy"}
 Toggle storage of invisible and escaped energy - by default this is false, and hence we do not store the invisible or escaped calibration hit truth (ctt) energy.
Gaudi::Property< bool > m_useCellWeights {this,"useCellWeights",false,"Toggle whether to use cell weights or not to calculate calibration hit contribution"}
 Toggle whether to use cell weights or not to calculate calibration hit contribution.

Detailed Description

Concrete tool that calculates calibration hit truth energies in xAOD::CaloCluster.

Definition at line 13 of file CaloCalibClusterTruthAttributerTool.h.

Constructor & Destructor Documentation

◆ CaloCalibClusterTruthAttributerTool()

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

Definition at line 7 of file CaloCalibClusterTruthAttributerTool.cxx.

7 : base_class(type,name,parent) {
8}

◆ ~CaloCalibClusterTruthAttributerTool()

CaloCalibClusterTruthAttributerTool::~CaloCalibClusterTruthAttributerTool ( )
virtualdefault

Member Function Documentation

◆ calculateTruthEnergies()

StatusCode CaloCalibClusterTruthAttributerTool::calculateTruthEnergies ( const xAOD::CaloCluster & theCaloCluster,
unsigned int numTruthParticles,
const std::map< Identifier, std::vector< const CaloCalibrationHit * > > & identifierToCaloHitMap,
std::vector< std::pair< unsigned int, double > > & barcodeTrueCalHitEnergy ) const
overridevirtual

This calculates the truth energies of N leading truth particles in a topocluster.

Definition at line 12 of file CaloCalibClusterTruthAttributerTool.cxx.

12 {
13
14 ATH_MSG_DEBUG("In calculateTruthEnergies");
15
16 const CaloClusterCellLink* theCellLinks = theCaloCluster.getCellLinks();
17
18 if (!theCellLinks) {
19 ATH_MSG_ERROR("A CaloCluster has no CaloClusterCellLinks");
20 return StatusCode::FAILURE;
21 }
22
23 std::map<unsigned int, double> truthIDTruePtMap;
24
25 //Loop on calorimeter cells to sum up the truth energies of the truth particles.
26 CaloClusterCellLink::const_iterator firstCell = theCellLinks->begin();
27 CaloClusterCellLink::const_iterator lastCell = theCellLinks->end();
28
29 for (; firstCell != lastCell; ++firstCell) {
30 const CaloCell* thisCaloCell = (*firstCell);
31
32 if (!thisCaloCell){
33 ATH_MSG_WARNING("Have invalid pointer to CaloCell");
34 continue;
35 }
36
37 //get the unique calorimeter cell identifier
38 Identifier cellID = thisCaloCell->ID();
39
40 //get the weight of the cell
41 double cellWeight = firstCell.weight();
42
43 //look up the calibration hit that corresponds to this calorimeter cell - we use find because not all calorimeter cells will have calibration hits
44 std::map<Identifier,std::vector<const CaloCalibrationHit*> >::const_iterator identifierToCaloHitMapIterator = identifierToCaloHitMap.find(cellID);
45 if (identifierToCaloHitMap.end() == identifierToCaloHitMapIterator) continue;
46 std::vector<const CaloCalibrationHit*> theseCalibrationHits = (*identifierToCaloHitMapIterator).second;
47
48 for (const auto *thisCalibrationHit : theseCalibrationHits){
49 const int truthID = HepMC::uniqueID(thisCalibrationHit);
50 double thisCalHitTruthEnergy = thisCalibrationHit->energyEM() + thisCalibrationHit->energyNonEM();
51 if (true == m_fullTruthEnergy) thisCalHitTruthEnergy += (thisCalibrationHit->energyEscaped() + thisCalibrationHit->energyInvisible());
52
53 //This only makes sense to use for clusters that are NOT calibrated
54 //If the cluster is calibrated this weight includes a calibration factor
55 //which is not relevant for truth information.
56 //For uncalibrated clusters it can contain a gemetrical weight
57 //and also a weight due to pflow reweighting of a cell energy
58 if (m_useCellWeights) thisCalHitTruthEnergy *= cellWeight;
59
60 if (truthID > HepMC::UNDEFINED_ID) {
61 truthIDTruePtMap[truthID] += thisCalHitTruthEnergy;
62 }
63
64 }//calibration hit loop
65
66 }//loop on calorimeter cells to sum up truth energies
67
68 //now create a vector with the same information as the map, which we can then sort
69 std::vector<std::pair<unsigned int, double > > truthIDTruePtPairs;
70
71 truthIDTruePtPairs.reserve(truthIDTruePtMap.size());
72 for (const auto& thisEntry : truthIDTruePtMap) truthIDTruePtPairs.emplace_back(thisEntry);
73
74 //sort vector by calibration hit truth energy
75 std::sort(truthIDTruePtPairs.begin(),truthIDTruePtPairs.end(),[]( std::pair<unsigned int, double> a, std::pair<unsigned int, double> b) -> bool {return a.second > b.second;} );
76
77 //store the truthID and truth energy of the top numTruthParticles truth particles
78 if (numTruthParticles > truthIDTruePtPairs.size()) numTruthParticles = truthIDTruePtPairs.size();
79 for ( unsigned int counter = 0; counter < numTruthParticles; counter++) truthIDTrueCalHitEnergy.push_back(truthIDTruePtPairs[counter]);
80
81 for (const auto& thisPair : truthIDTrueCalHitEnergy) ATH_MSG_DEBUG("Truncated loop 2: truthID and true energy are " << thisPair.first << " and " << thisPair.second << " for cluster with e, eta of " << theCaloCluster.e() << " and " << theCaloCluster.eta() );
82
83 return StatusCode::SUCCESS;
84
85}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
static Double_t a
Gaudi::Property< bool > m_useCellWeights
Toggle whether to use cell weights or not to calculate calibration hit contribution.
Gaudi::Property< bool > m_fullTruthEnergy
Toggle storage of invisible and escaped energy - by default this is false, and hence we do not store ...
Identifier ID() const
get ID (from cached data member) non-virtual and inline for fast access
Definition CaloCell.h:295
const CaloClusterCellLink * getCellLinks() const
Get a pointer to the CaloClusterCellLink object (const version)
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double e() const
The total energy of the particle.
int uniqueID(const T &p)
constexpr int UNDEFINED_ID
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.

Member Data Documentation

◆ m_fullTruthEnergy

Gaudi::Property<bool> CaloCalibClusterTruthAttributerTool::m_fullTruthEnergy {this,"storeFullTruthEnergy",false,"Toggle storage of invisible and escaped energy"}
private

Toggle storage of invisible and escaped energy - by default this is false, and hence we do not store the invisible or escaped calibration hit truth (ctt) energy.

Hence only the EM and nonEM truth ctt energy is stored by default

Definition at line 28 of file CaloCalibClusterTruthAttributerTool.h.

28{this,"storeFullTruthEnergy",false,"Toggle storage of invisible and escaped energy"};

◆ m_useCellWeights

Gaudi::Property<bool> CaloCalibClusterTruthAttributerTool::m_useCellWeights {this,"useCellWeights",false,"Toggle whether to use cell weights or not to calculate calibration hit contribution"}
private

Toggle whether to use cell weights or not to calculate calibration hit contribution.

Definition at line 31 of file CaloCalibClusterTruthAttributerTool.h.

31{this,"useCellWeights",false,"Toggle whether to use cell weights or not to calculate calibration hit contribution"};

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