ATLAS Offline Software
Loading...
Searching...
No Matches
CalibHitToCaloCellTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7
8// Calo include
14#include "CaloDetDescr/CaloDetDescrElement.h"
15
21
22#include "LArRecEvent/LArCell.h"
23#include "TileEvent/TileCell.h"
25
26CalibHitToCaloCellTool::CalibHitToCaloCellTool(const std::string& t, const std::string& n, const IInterface* p)
27 : AthAlgTool(t,n,p),
28 m_caloGain((int)CaloGain::LARLOWGAIN),
29 m_caloCell_Tot("TotalCalibCell"), m_caloCell_Vis("VisCalibCell"),
31{
32 declareInterface<CalibHitToCaloCellTool>(this);
33
34 declareProperty("CaloGain", m_caloGain);
35 declareProperty("CalibHitContainers", m_calibHitContainerNames);
36
37 declareProperty("CellTotEne", m_caloCell_Tot);
38 declareProperty("CellVisEne", m_caloCell_Vis);
39 declareProperty("CellEmEne", m_caloCell_Em);
40 declareProperty("CellNonEmEne", m_caloCell_NonEm);
41 declareProperty("DoTile", m_doTile=false);
42
43 declareProperty("OutputCellContainerName", m_outputCellContainerName = "TruthCells");
44 declareProperty("OutputClusterContainerName", m_outputClusterContainerName = "TruthClusters");
45
46 m_tileActiveHitCnt = "TileCalibHitActiveCell";
47 m_tileInactiveHitCnt = "TileCalibHitInactiveCell";
48 m_tileDMHitCnt = "TileCalibHitDeadMaterial";
49 m_larActHitCnt = "LArCalibrationHitActive";
50 m_larInactHitCnt = "LArCalibrationHitInactive";
51 m_larDMHitCnt = "LArCalibrationHitDeadMaterial";
52
53}
54
55
57
58
61{
62 // retrieve ID helpers from det store
63 ATH_MSG_INFO("initialisation ID helpers" );
64
65 ATH_CHECK( detStore()->retrieve(m_caloCell_ID) );
66 ATH_CHECK( detStore()->retrieve(m_caloDM_ID) );
67 ATH_CHECK( m_caloMgrKey.initialize() );
68
69
70 for (unsigned int i=0; i<CalibHitUtils::nEnergyTypes; i++) {
74 }
75
76 ATH_CHECK(m_cellContKeys.initialize());
77 ATH_CHECK(m_clusterContKeys.initialize());
78 ATH_CHECK(m_cellLinkKeys.initialize());
79
80 ATH_MSG_INFO("initialisation completed" );
81 return StatusCode::SUCCESS;
82}
83
84
87{
88 ATH_MSG_DEBUG("in calibHitToCaloCellTool");
89
91 ATH_CHECK(caloMgrHandle.isValid());
92 const CaloDetDescrManager* caloDDMgr = *caloMgrHandle;
93
94 //CaloCellContainer* truthCells[3];
95 //xAOD::CaloClusterContainer* truthClusters[3];
96
97 std::vector<SG::WriteHandle<CaloCellContainer> > truthCells;
98 std::vector<SG::WriteHandle<xAOD::CaloClusterContainer> > truthClusters;
99 std::vector<SG::WriteHandle<CaloClusterCellLinkContainer> > truthLinks;
100
101
102 // register containers for cells and clusters
103 for (unsigned int i=0; i<CalibHitUtils::nEnergyTypes; i++) {
104
105 truthCells.emplace_back(m_cellContKeys[i]);
106 ATH_CHECK(truthCells.back().record(std::make_unique<CaloCellContainer>()));
107
108 truthClusters.emplace_back(m_clusterContKeys[i]);
110
111 truthLinks.emplace_back(m_cellLinkKeys[i]);
112 ATH_CHECK(truthLinks.back().record(std::make_unique<CaloClusterCellLinkContainer>()));
113 }
114
115 // retrieve calibration hit containers
116 const unsigned int nCont = m_calibHitContainerNames.size();
117 std::vector<const CaloCalibrationHitContainer*> calibHitContainers(nCont,nullptr);
118 for (unsigned int i=0; i<nCont; i++) {
119 ATH_MSG_DEBUG("Retrieving " << m_calibHitContainerNames[i]);
120 ATH_CHECK( evtStore()->retrieve(calibHitContainers[i], m_calibHitContainerNames[i].c_str()) );
121 ATH_MSG_DEBUG(" Retrieved container " << calibHitContainers[i]->Name() << " with size " << calibHitContainers[i]->Size() );
122 if( calibHitContainers[i]->Size() == 0 ) {
123 ATH_MSG_DEBUG("Container " << calibHitContainers[i]->Name() << " is empty");
124 }
125 }
126 ATH_MSG_DEBUG("CaloCalibrationHitContainers retrieved successfuly" );
127
128 //count
129 int nchan=0;
130 int em_nchan=0;
131 int hec_nchan=0;
132 int fcal_nchan=0;
133 int tile_nchan=0;
134
135 std::vector<Identifier> ID;
136
137
138 std::vector<CaloCell*> CellsEtot;
139 std::vector<CaloCell*> CellsEvis;
140 std::vector<CaloCell*> CellsEem;
141
142 int nhitsInactive = 0;
143
144 for (unsigned int i=0; i<calibHitContainers.size(); i++) {
145 for( const auto *const calibhit: *(calibHitContainers[i])) {
146 //care only for deposits of the given truth particle
147 if (!MC::isGenStable(HepMC::barcode(calibhit))) continue; // FIXME barcode-based
148
149 double Etot = calibhit->energyTotal();
150 double Eem = calibhit->energy(0);
151 double Enonem = calibhit->energy(1);
152 double Evis = Eem + Enonem;
153
154 Identifier id=calibhit->cellID();
155
156 //merge inactive and active hits from the same cell together
157 if (i>0) {
158 //find if have already created cell..
159 bool isNewId = true;
160 for (int n=0; n<nhitsInactive; n++) {
161 if( id == ID[n] ) { //found
162 CellsEtot[n]->addEnergy(Etot);
163 CellsEvis[n]->addEnergy(Evis);
164 CellsEem[n]->addEnergy(Eem);
165 isNewId = false;
166 break;
167 }
168 }
169 if(!isNewId) continue; //go to next hit, else create new cell for this hit
170 }
171
172 //check if this ID is LAr or Tile
173 if(m_caloCell_ID->is_lar(id)) {
174 ATH_MSG_VERBOSE( "Found LAr cell" );
175 const CaloDetDescrElement* caloDDE = caloDDMgr->get_element(id);
176 CellsEtot.push_back(new LArCell(caloDDE, id, Etot, 0., 0, 0, (CaloGain::CaloGain)m_caloGain)) ;
177 CellsEvis.push_back(new LArCell(caloDDE, id, Evis, 0., 0, 0, (CaloGain::CaloGain)m_caloGain));
178 CellsEem.push_back(new LArCell(caloDDE, id, Eem, 0., 0, 0, (CaloGain::CaloGain)m_caloGain));
179 ID.push_back(id);
180 ++nchan;
181 }
182 else if(m_caloCell_ID->is_tile(id)) {
183 ATH_MSG_VERBOSE( "Found Tile cell" );
184 const CaloDetDescrElement* caloDDE = caloDDMgr->get_element(id);
185 CellsEtot.push_back(new TileCell(caloDDE, id, Etot, 0., 0, 0, (CaloGain::CaloGain)m_caloGain)) ;
186 CellsEvis.push_back(new TileCell(caloDDE, id, Evis, 0., 0, 0, (CaloGain::CaloGain)m_caloGain));
187 CellsEem.push_back(new TileCell(caloDDE, id, Eem, 0., 0, 0, (CaloGain::CaloGain)m_caloGain));
188 ID.push_back(id);
189 ++nchan;
190 }
191 else { //other, DeadMaterial
193 ATH_MSG_VERBOSE( "Found unknown cell" );
194 continue;
195 }
196 }
197 if (i==0) nhitsInactive = (int)ID.size();
198 }
199
200 //Now, put cells in the containers keeping the order. First goes EM, then HEC and so on
201 // if(CellsEtot.size()==0) {
202 // ID.clear();
203 // return StatusCode::SUCCESS;
204 // }
205
206 ATH_MSG_DEBUG("N cells : " << nchan );
207
208 for(int itr=0; itr!=nchan; itr++) {
209 if(m_caloCell_ID->is_em(CellsEtot[itr]->ID())) {
210 truthCells[CalibHitUtils::EnergyTotal]->push_back(CellsEtot[itr]);
211 truthCells[CalibHitUtils::EnergyVisible]->push_back(CellsEvis[itr]);
212 truthCells[CalibHitUtils::EnergyEM]->push_back(CellsEem[itr]);
213 ++em_nchan;
214 }
215 }
216 if(em_nchan) {
217 for (int i=0;i<CalibHitUtils::nEnergyTypes;i++) truthCells[i]->setHasCalo(CaloCell_ID::LAREM);
218 }
219
220 for(int itr=0; itr!=nchan; itr++) {
221 if(m_caloCell_ID->is_hec(CellsEtot[itr]->ID())) {
222 truthCells[CalibHitUtils::EnergyTotal]->push_back(CellsEtot[itr]);
223 truthCells[CalibHitUtils::EnergyVisible]->push_back(CellsEvis[itr]);
224 truthCells[CalibHitUtils::EnergyEM]->push_back(CellsEem[itr]);
225 ++hec_nchan;
226 }
227 }
228 if(hec_nchan){
229 for (int i=0;i<CalibHitUtils::nEnergyTypes;i++) truthCells[i]->setHasCalo(CaloCell_ID::LARHEC);
230 }
231
232 for(int itr=0; itr!=nchan; itr++) {
233 if(m_caloCell_ID->is_fcal(CellsEtot[itr]->ID())) {
234 truthCells[CalibHitUtils::EnergyTotal]->push_back(CellsEtot[itr]);
235 truthCells[CalibHitUtils::EnergyVisible]->push_back(CellsEvis[itr]);
236 truthCells[CalibHitUtils::EnergyEM]->push_back(CellsEem[itr]);
237 ++fcal_nchan;
238 }
239 }
240 if(fcal_nchan) {
241 for (int i=0;i<CalibHitUtils::nEnergyTypes;i++) truthCells[i]->setHasCalo(CaloCell_ID::LARFCAL);
242 }
243
244 for(int itr=0; itr!=nchan; itr++) {
245 if((m_caloCell_ID->is_tile(CellsEtot[itr]->ID()))) {
246 truthCells[CalibHitUtils::EnergyTotal]->push_back(CellsEtot[itr]);
247 truthCells[CalibHitUtils::EnergyVisible]->push_back(CellsEvis[itr]);
248 truthCells[CalibHitUtils::EnergyEM]->push_back(CellsEem[itr]);
249 ++tile_nchan;
250 }
251 }
252 if(tile_nchan) {
253 for (int i=0;i<CalibHitUtils::nEnergyTypes;i++) truthCells[i]->setHasCalo(CaloCell_ID::TILE);
254 }
255 ATH_MSG_DEBUG("--- LAr INFO --- "<<nchan );
256 ATH_MSG_DEBUG("LArCells = "<<nchan );
257 ATH_MSG_DEBUG("EMCells = "<<em_nchan );
258 ATH_MSG_DEBUG("HECCells = "<<hec_nchan );
259 ATH_MSG_DEBUG("FCALCells = "<<fcal_nchan );
260 ATH_MSG_DEBUG("TileCells = "<<tile_nchan );
261
262 ID.clear();
263
265
266 ATH_MSG_DEBUG("making truth cluster");
267 xAOD::CaloCluster* truthCluster[3] = {nullptr,nullptr,nullptr};
268 for (int i=0;i<CalibHitUtils::nEnergyTypes;i++) {
269 truthCluster[i] = CaloClusterStoreHelper::makeCluster(truthClusters[i].ptr(),truthCells[i].ptr());
270 if (!truthCluster[i]) {
271 ATH_MSG_FATAL("makeCluster failed");
272 return StatusCode::FAILURE;
273 }
274 for (const auto cell: *truthCells[i]) {
275 if(m_caloCell_ID->is_lar(cell->ID()) || m_caloCell_ID->is_tile(cell->ID()))
276 truthCluster[i]->addCell( truthCells[i]->findIndex(cell->caloDDE()->calo_hash()) , 1.);
277 }
278
280 CaloClusterKineHelper::calculateKine(truthCluster[i], true, true);
282 truthClusters[i].ptr()));
283
284 ATH_MSG_INFO("Created truth cluster with " << m_energyTypeToStr[i] <<" " << truthCluster[i]->e());
285 }
286
287 ATH_MSG_DEBUG("execute() completed successfully" );
288 return StatusCode::SUCCESS;
289}
290
291
294{
295 ATH_MSG_INFO("finalize() successfully" );
296 return StatusCode::SUCCESS;
297}
298
#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)
std::vector< Identifier > ID
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:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
ServiceHandle< StoreGateSvc > & evtStore()
const ServiceHandle< StoreGateSvc > & detStore() const
SG::WriteHandleKeyArray< CaloCellContainer > m_cellContKeys
SG::WriteHandleKeyArray< CaloClusterCellLinkContainer > m_cellLinkKeys
SG::ReadCondHandleKey< CaloDetDescrManager > m_caloMgrKey
CalibHitToCaloCellTool(const std::string &t, const std::string &n, const IInterface *p)
const std::array< std::string, 3 > m_energyTypeToStr
std::vector< std::string > m_calibHitContainerNames
SG::WriteHandleKeyArray< xAOD::CaloClusterContainer > m_clusterContKeys
StatusCode processCalibHitsFromParticle() const
const CaloCell_ID * m_caloCell_ID
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.