ATLAS Offline Software
Loading...
Searching...
No Matches
CalibHitToCaloCellTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7// Calo include
13#include "CaloDetDescr/CaloDetDescrElement.h"
14
20
21#include "LArRecEvent/LArCell.h"
22#include "TileEvent/TileCell.h"
24
25#include <memory>
26
27CalibHitToCaloCellTool::CalibHitToCaloCellTool(const std::string& t, const std::string& n, const IInterface* p)
28 : AthAlgTool(t,n,p)
29{
30 declareInterface<CalibHitToCaloCellTool>(this);
31}
32
34
37{
38 // retrieve ID helpers from det store
39 ATH_MSG_INFO("initialisation ID helpers" );
40
41 ATH_CHECK( detStore()->retrieve(m_caloCell_ID) );
42 ATH_CHECK( detStore()->retrieve(m_caloDM_ID) );
43 ATH_CHECK( m_caloMgrKey.initialize() );
44
45
46 for (unsigned int i=0; i<CalibHitUtils::nEnergyTypes; i++) {
50 }
51
52 ATH_CHECK(m_cellContKeys.initialize());
53 ATH_CHECK(m_clusterContKeys.initialize());
54 ATH_CHECK(m_cellLinkKeys.initialize());
55
56 ATH_MSG_INFO("initialisation completed" );
57 return StatusCode::SUCCESS;
58}
59
60
62{
63 ATH_MSG_DEBUG("in calibHitToCaloCellTool");
64
66 ATH_CHECK(caloMgrHandle.isValid());
67 const CaloDetDescrManager* caloDDMgr = *caloMgrHandle;
68
69 std::vector<SG::WriteHandle<CaloCellContainer> > truthCells;
70 std::vector<SG::WriteHandle<xAOD::CaloClusterContainer> > truthClusters;
71 std::vector<SG::WriteHandle<CaloClusterCellLinkContainer> > truthLinks;
72
73
74 // register containers for cells and clusters
75 for (unsigned int i=0; i<CalibHitUtils::nEnergyTypes; i++) {
76
77 truthCells.emplace_back(m_cellContKeys[i]);
78 ATH_CHECK(truthCells.back().record(std::make_unique<CaloCellContainer>()));
79
80 truthClusters.emplace_back(m_clusterContKeys[i]);
82
83 truthLinks.emplace_back(m_cellLinkKeys[i]);
84 ATH_CHECK(truthLinks.back().record(std::make_unique<CaloClusterCellLinkContainer>()));
85 }
86
87 // retrieve calibration hit containers
88 const unsigned int nCont = m_calibHitContainerNames.size();
89 std::vector<const CaloCalibrationHitContainer*> calibHitContainers(nCont,nullptr);
90 for (unsigned int i=0; i<nCont; i++) {
91 ATH_MSG_DEBUG("Retrieving " << m_calibHitContainerNames[i]);
92 ATH_CHECK( evtStore()->retrieve(calibHitContainers[i], m_calibHitContainerNames[i].c_str()) );
93 ATH_MSG_DEBUG(" Retrieved container " << calibHitContainers[i]->Name() << " with size " << calibHitContainers[i]->Size() );
94 if( calibHitContainers[i]->Size() == 0 ) {
95 ATH_MSG_DEBUG("Container " << calibHitContainers[i]->Name() << " is empty");
96 }
97 }
98 ATH_MSG_DEBUG("CaloCalibrationHitContainers retrieved successfuly" );
99
100 //count
101 int nchan=0;
102 int em_nchan=0;
103 int hec_nchan=0;
104 int fcal_nchan=0;
105 int tile_nchan=0;
106
107 std::vector<Identifier> ID;
108
109 using CellPtr = std::unique_ptr<CaloCell>;
110 std::vector<CellPtr> CellsEtot;
111 std::vector<CellPtr> CellsEvis;
112 std::vector<CellPtr> CellsEem;
113
114 int nhitsInactive = 0;
115
116 for (unsigned int i=0; i<calibHitContainers.size(); i++) {
117 for( const auto *const calibhit: *(calibHitContainers[i])) {
118 //care only for deposits of the given truth particle
119 if (!MC::isGenStable(HepMC::barcode(calibhit))) continue; // FIXME barcode-based
120
121 double Etot = calibhit->energyTotal();
122 double Eem = calibhit->energy(0);
123 double Enonem = calibhit->energy(1);
124 double Evis = Eem + Enonem;
125
126 Identifier id=calibhit->cellID();
127
128 //merge inactive and active hits from the same cell together
129 if (i>0) {
130 //find if have already created cell..
131 bool isNewId = true;
132 for (int n=0; n<nhitsInactive; n++) {
133 if( id == ID[n] ) { //found
134 CellsEtot[n]->addEnergy(Etot);
135 CellsEvis[n]->addEnergy(Evis);
136 CellsEem[n]->addEnergy(Eem);
137 isNewId = false;
138 break;
139 }
140 }
141 if(!isNewId) continue; //go to next hit, else create new cell for this hit
142 }
143
144 //check if this ID is LAr or Tile
145 if(m_caloCell_ID->is_lar(id)) {
146 ATH_MSG_VERBOSE( "Found LAr cell" );
147 const CaloDetDescrElement* caloDDE = caloDDMgr->get_element(id);
148 CellsEtot.push_back(CellPtr{new LArCell(caloDDE, id, Etot, 0., 0, 0, (CaloGain::CaloGain)m_caloGain.value())});
149 CellsEvis.push_back(CellPtr{new LArCell(caloDDE, id, Evis, 0., 0, 0, (CaloGain::CaloGain)m_caloGain.value())});
150 CellsEem.push_back(CellPtr{new LArCell(caloDDE, id, Eem, 0., 0, 0, (CaloGain::CaloGain)m_caloGain.value())});
151 ID.push_back(id);
152 ++nchan;
153 }
154 else if(m_caloCell_ID->is_tile(id)) {
155 ATH_MSG_VERBOSE( "Found Tile cell" );
156 const CaloDetDescrElement* caloDDE = caloDDMgr->get_element(id);
157 CellsEtot.push_back(CellPtr{new TileCell(caloDDE, id, Etot, 0., 0, 0, (CaloGain::CaloGain)m_caloGain.value())});
158 CellsEvis.push_back(CellPtr{new TileCell(caloDDE, id, Evis, 0., 0, 0, (CaloGain::CaloGain)m_caloGain.value())});
159 CellsEem.push_back(CellPtr{new TileCell(caloDDE, id, Eem, 0., 0, 0, (CaloGain::CaloGain)m_caloGain.value())});
160 ID.push_back(id);
161 ++nchan;
162 }
163 else { //other, DeadMaterial
165 ATH_MSG_VERBOSE( "Found unknown cell" );
166 continue;
167 }
168 }
169 if (i==0) nhitsInactive = (int)ID.size();
170 }
171
172 //Now, put cells in the containers keeping the order. First goes EM, then HEC and so on
173 // if(CellsEtot.size()==0) {
174 // ID.clear();
175 // return StatusCode::SUCCESS;
176 // }
177
178 ATH_MSG_DEBUG("N cells : " << nchan );
179
180 for(int itr=0; itr!=nchan; itr++) {
181 if(m_caloCell_ID->is_em(CellsEtot[itr]->ID())) {
182 truthCells[CalibHitUtils::EnergyTotal]->push_back(std::move(CellsEtot[itr]));
183 truthCells[CalibHitUtils::EnergyVisible]->push_back(std::move(CellsEvis[itr]));
184 truthCells[CalibHitUtils::EnergyEM]->push_back(std::move(CellsEem[itr]));
185 ++em_nchan;
186 }
187 }
188 if(em_nchan) {
189 for (int i=0;i<CalibHitUtils::nEnergyTypes;i++) truthCells[i]->setHasCalo(CaloCell_ID::LAREM);
190 }
191
192 for(int itr=0; itr!=nchan; itr++) {
193 if(m_caloCell_ID->is_hec(CellsEtot[itr]->ID())) {
194 truthCells[CalibHitUtils::EnergyTotal]->push_back(std::move(CellsEtot[itr]));
195 truthCells[CalibHitUtils::EnergyVisible]->push_back(std::move(CellsEvis[itr]));
196 truthCells[CalibHitUtils::EnergyEM]->push_back(std::move(CellsEem[itr]));
197 ++hec_nchan;
198 }
199 }
200 if(hec_nchan){
201 for (int i=0;i<CalibHitUtils::nEnergyTypes;i++) truthCells[i]->setHasCalo(CaloCell_ID::LARHEC);
202 }
203
204 for(int itr=0; itr!=nchan; itr++) {
205 if(m_caloCell_ID->is_fcal(CellsEtot[itr]->ID())) {
206 truthCells[CalibHitUtils::EnergyTotal]->push_back(std::move(CellsEtot[itr]));
207 truthCells[CalibHitUtils::EnergyVisible]->push_back(std::move(CellsEvis[itr]));
208 truthCells[CalibHitUtils::EnergyEM]->push_back(std::move(CellsEem[itr]));
209 ++fcal_nchan;
210 }
211 }
212 if(fcal_nchan) {
213 for (int i=0;i<CalibHitUtils::nEnergyTypes;i++) truthCells[i]->setHasCalo(CaloCell_ID::LARFCAL);
214 }
215
216 for(int itr=0; itr!=nchan; itr++) {
217 if((m_caloCell_ID->is_tile(CellsEtot[itr]->ID()))) {
218 truthCells[CalibHitUtils::EnergyTotal]->push_back(std::move(CellsEtot[itr]));
219 truthCells[CalibHitUtils::EnergyVisible]->push_back(std::move(CellsEvis[itr]));
220 truthCells[CalibHitUtils::EnergyEM]->push_back(std::move(CellsEem[itr]));
221 ++tile_nchan;
222 }
223 }
224 if(tile_nchan) {
225 for (int i=0;i<CalibHitUtils::nEnergyTypes;i++) truthCells[i]->setHasCalo(CaloCell_ID::TILE);
226 }
227 ATH_MSG_DEBUG("--- LAr INFO --- "<<nchan );
228 ATH_MSG_DEBUG("LArCells = "<<nchan );
229 ATH_MSG_DEBUG("EMCells = "<<em_nchan );
230 ATH_MSG_DEBUG("HECCells = "<<hec_nchan );
231 ATH_MSG_DEBUG("FCALCells = "<<fcal_nchan );
232 ATH_MSG_DEBUG("TileCells = "<<tile_nchan );
233
234 ID.clear();
235
237
238 ATH_MSG_DEBUG("making truth cluster");
239 xAOD::CaloCluster* truthCluster[3] = {nullptr,nullptr,nullptr};
240 for (int i=0;i<CalibHitUtils::nEnergyTypes;i++) {
241 truthCluster[i] = CaloClusterStoreHelper::makeCluster(truthClusters[i].ptr(),truthCells[i].ptr());
242 if (!truthCluster[i]) {
243 ATH_MSG_FATAL("makeCluster failed");
244 return StatusCode::FAILURE;
245 }
246 for (const auto cell: *truthCells[i]) {
247 if(m_caloCell_ID->is_lar(cell->ID()) || m_caloCell_ID->is_tile(cell->ID()))
248 truthCluster[i]->addCell( truthCells[i]->findIndex(cell->caloDDE()->calo_hash()) , 1.);
249 }
250
252 CaloClusterKineHelper::calculateKine(truthCluster[i], true, true);
254 truthClusters[i].ptr()));
255
256 ATH_MSG_INFO("Created truth cluster with " << m_energyTypeToStr[i] <<" " << truthCluster[i]->e());
257 }
258
259 ATH_MSG_DEBUG("execute() completed successfully" );
260 return StatusCode::SUCCESS;
261}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_DEBUG(x)
Convert energy deposits from calibration hits to CaloCell, xAOD::CaloCluster.
ATLAS-specific HepMC functions.
JetDumper::Name Name
Definition JetDumper.cxx:19
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
ServiceHandle< StoreGateSvc > & evtStore()
const ServiceHandle< StoreGateSvc > & detStore() const
SG::WriteHandleKeyArray< CaloCellContainer > m_cellContKeys
SG::WriteHandleKeyArray< CaloClusterCellLinkContainer > m_cellLinkKeys
Gaudi::Property< std::string > m_outputCellContainerName
SG::ReadCondHandleKey< CaloDetDescrManager > m_caloMgrKey
Gaudi::Property< std::vector< std::string > > m_calibHitContainerNames
Gaudi::Property< int > m_caloGain
CalibHitToCaloCellTool(const std::string &t, const std::string &n, const IInterface *p)
const std::array< std::string, 3 > m_energyTypeToStr
SG::WriteHandleKeyArray< xAOD::CaloClusterContainer > m_clusterContKeys
Gaudi::Property< std::string > m_outputClusterContainerName
StatusCode processCalibHitsFromParticle() const
const CaloCell_ID * m_caloCell_ID
virtual StatusCode initialize() override
static void calculateKine(xAOD::CaloCluster *clu, const bool useweight=true, const bool updateLayers=true, const bool useGPUCriteria=false)
Helper class to calculate cluster kinematics based on cells.
static std::unique_ptr< xAOD::CaloCluster > makeCluster(const CaloCellContainer *cellCont)
Creates a valid CaloCluster with a private Aux-Store and CellLink container.
static StatusCode AddContainerWriteHandle(SG::WriteHandle< xAOD::CaloClusterContainer > &clusColl)
Creates a new xAOD::CaloClusterContainer in the given WriteHandle + CaloClusterAuxContainer and recor...
static StatusCode finalizeClusters(SG::WriteHandle< CaloClusterCellLinkContainer > &h, xAOD::CaloClusterContainer *pClusterColl)
Finalize clusters (move CaloClusterCellLink to a separate container).
This class groups all DetDescr information related to a CaloCell.
const CaloDetDescrElement * get_element(const Identifier &cellId) const
get element by its identifier
This class provides the client interface for accessing the detector description information common to...
Data object for LAr calorimeter readout cell.
Definition LArCell.h:53
void setClusterSize(const ClusterSize)
Get cluster size.
bool addCell(const unsigned index, const double weight)
Method to add a cell to the cluster (Beware: Kinematics not updated!).
int barcode(const T *p)
Definition Barcode.h:16
bool isGenStable(const T &p)
Determine if the particle is stable at the generator (not det-sim) level,.
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.