ATLAS Offline Software
Loading...
Searching...
No Matches
LArFCS_StepInfoSD.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
5// class header
6#include "LArFCS_StepInfoSD.h"
8#include "G4Step.hh"
9#include "G4ThreeVector.hh"
10#include "G4TouchableHistory.hh"
11
13#include <string>
14#include <utility>
15
16#include "CaloDetDescr/CaloDetDescrElement.h"
18#include "GaudiKernel/MsgStream.h"
20
23 : FCS_StepInfoSD(std::move(a_name), config),
24 m_calculator(m_config.m_LArCalculator) {}
25
27
28G4bool LArFCS_StepInfoSD::ProcessHits(G4Step* a_step, G4TouchableHistory*) {
29 G4bool result(false);
30 // If there's no energy, there's no hit. (Aside: Isn't this energy
31 // the same as the energy from the calculator? Not necessarily.
32 // The calculator may include detector effects such as
33 // charge-collection which are not modeled by Geant4.)
34 if (a_step->GetTotalEnergyDeposit() <= 0.) {
35 return result;
36 }
37
38 if (m_calculator) {
39 if (!m_calo_dd_man.get()) {
41 }
42
43 const double StepLength = a_step->GetStepLength() / CLHEP::mm;
44 const G4ThreeVector preStepPosition =
45 a_step->GetPreStepPoint()
46 ->GetPosition(); // pre step is the position we're interested in
47 const G4ThreeVector postStepPosition =
48 a_step->GetPostStepPoint()->GetPosition();
49 std::vector<const G4Step*> steps;
50 bool madeSubSteps(false);
51 if (m_config.shorten_lar_step && StepLength > m_config.substpsize) {
52 // create smaller substeps instead
53 G4int nsub_step = (int)(StepLength / m_config.substpsize) + 1;
54 G4double delta = 1. / ((double)nsub_step);
55 // G4cout <<"Orig prestep: "<<preStepPosition<<std::endl;
56 for (G4int i = 0; i < nsub_step; i++) {
57 G4double fraction1 = ((G4double)i) * delta;
58 G4double fraction2 = (((G4double)i) + 1.) * delta;
59 G4ThreeVector subpoint1 =
60 preStepPosition * (1 - fraction1) + postStepPosition * fraction1;
61 G4ThreeVector subpoint2 =
62 preStepPosition * (1 - fraction2) + postStepPosition * fraction2;
63 G4StepPoint* startpoint = new G4StepPoint(*(a_step->GetPreStepPoint()));
64 G4StepPoint* endpoint = new G4StepPoint(*(a_step->GetPostStepPoint()));
65 startpoint->SetPosition(subpoint1);
66 endpoint->SetPosition(subpoint2);
67
68 G4Step* newstep = new G4Step(*a_step);
69 newstep->SetPreStepPoint(startpoint);
70 newstep->SetPostStepPoint(endpoint);
71 newstep->SetStepLength((subpoint1 - subpoint2).mag());
72 newstep->SetTotalEnergyDeposit(a_step->GetTotalEnergyDeposit() /
73 nsub_step);
74 steps.push_back(newstep);
75 madeSubSteps = true;
76 }
77 } else {
78 steps.push_back(a_step);
79 }
80
81 if (m_config.verboseLevel > 4) {
82 G4cout << "LArFCS_StepInfoSD::ProcessHits(): original step in Volume: "
83 << a_step->GetPreStepPoint()->GetPhysicalVolume()->GetName()
84 << " position: " << a_step->GetPreStepPoint()->GetPosition()
85 << " Length=" << a_step->GetStepLength() / CLHEP::mm
86 << " E=" << a_step->GetTotalEnergyDeposit() << G4endl;
87 std::vector<LArHitData> processedHits;
88 if (m_calculator->Process(a_step, processedHits)) {
89 if (m_config.verboseLevel > -1) {
90 G4cout << " #LArHitData=" << processedHits.size();
91 for (const auto& lhd : processedHits) {
92 G4cout << " ; id=" << (std::string)lhd.id << " E=" << lhd.energy
93 << G4endl;
94 }
95 G4cout << G4endl;
96 }
97 }
98 G4cout << "LArFCS_StepInfoSD::ProcessHits(): #substep=" << steps.size()
99 << G4endl;
100 }
101
102 double et_all = 0; // Total collected charge in all substeps
103 for (const G4Step* substep : steps) {
104 double et(0.); // Total collected charge in this substep
105 G4ThreeVector stepPosition =
106 0.5 * (substep->GetPreStepPoint()->GetPosition() +
107 substep->GetPostStepPoint()->GetPosition());
108 std::vector<LArHitData> processedHits;
109 if (m_config.verboseLevel > 4) {
110 G4cout << " LArFCS_StepInfoSD::ProcessHits(): substep in Volume: "
111 << substep->GetPreStepPoint()->GetPhysicalVolume()->GetName()
112 << " position: " << substep->GetPreStepPoint()->GetPosition()
113 << " Length=" << substep->GetStepLength() / CLHEP::mm
114 << " E=" << substep->GetTotalEnergyDeposit() << G4endl;
115 }
116 if (m_calculator->Process(substep, processedHits)) {
117 for (const auto& larhit : processedHits) {
118 et += larhit.energy;
119 et_all += larhit.energy;
120 }
121 if (m_config.verboseLevel > 4) {
122 G4cout << " substep #LArHitData=" << processedHits.size();
123 for (const auto& lhd : processedHits) {
124 G4cout << " ; id=" << (std::string)lhd.id << " E=" << lhd.energy
125 << G4endl;
126 }
127 G4cout << G4endl;
128 }
129 } else {
130 if (m_config.verboseLevel > 4) {
131 // Maybe 0 hits or something like that...
132 G4cout << this->GetName()
133 << " WARNING ProcessHits: Call to ILArCalculatorSvc::Process "
134 "failed! Details:"
135 << G4endl << " " << "Volume: "
136 << a_step->GetPreStepPoint()->GetPhysicalVolume()->GetName()
137 << " " << m_calculator << " position: " << stepPosition
138 << " SL: " << StepLength << G4endl << " "
139 << "Orig position: "
140 << substep->GetPreStepPoint()->GetPosition() << " / "
141 << substep->GetPostStepPoint()->GetPosition() << G4endl
142 << " " << "StepLength: " << StepLength
143 << " step: " << preStepPosition << " / " << postStepPosition
144 << G4endl;
145 }
146 continue;
147 }
148
149 // drop hits with zero deposited energy (could still happen with negative
150 // corrections from calculator)
151 // Or if total energy is <0
152 if (et <= 0.) {
153 if (m_config.verboseLevel > 4) {
154 G4cout << this->GetName()
155 << " WARNING ProcessHits: Total negative energy: " << et
156 << " not processing..." << G4endl;
157 }
158 continue;
159 }
160
161 const size_t numberOfProcessedHits = processedHits.size();
162 const G4ThreeVector originalStepPosition = stepPosition;
163 double maxSubHitEnergy = -999.;
164 int maxSubHitEnergyindex = -1;
165 if (numberOfProcessedHits > 0) {
166 maxSubHitEnergy = processedHits[0].energy;
167 maxSubHitEnergyindex = 0;
168 }
169 // Figure out the subhit with most energy
170 for (size_t i = 1; i < numberOfProcessedHits; ++i) {
171 if (maxSubHitEnergy < processedHits[i].energy) {
172 maxSubHitEnergy = processedHits[i].energy;
173 maxSubHitEnergyindex = i;
174 }
175 }
176 if (maxSubHitEnergyindex ==
177 -1) { // because there were no hits; numberOfProcessedHits ==0
178 G4cout << this->GetName()
179 << " WARNING ProcessHits: numberOfProcessedHits is zero"
180 << G4endl;
181 continue;
182 }
183 // Identifier for the subhit with max energy
184 Identifier maxEnergyIdentifier =
185 this->ConvertID(processedHits[maxSubHitEnergyindex].id);
186 const CaloDetDescrElement* maxEnergyCell =
187 m_calo_dd_man.get()->get_element(maxEnergyIdentifier);
188
189 Identifier invalidIdentifier;
190 // for (size_t i=0; i<numberOfProcessedHits; ++i) {
191 for (const auto& larhit : processedHits) {
192 Identifier id = this->ConvertID(larhit.id);
193 if (larhit.id[0] == 10) {
194 if (m_config.verboseLevel > 9) {
195 G4cout << this->GetName()
196 << " VERBOSE ProcessHits: Dead Material LArG4Identifier: "
197 << (std::string)larhit.id << G4endl << " " << id
198 << G4endl << " " << id.getString() << G4endl
199 << " "
200 << a_step->GetPreStepPoint()
201 ->GetPhysicalVolume()
202 ->GetLogicalVolume()
203 ->GetName()
204 << G4endl << " numberOfProcessedHits: "
205 << numberOfProcessedHits << G4endl;
206 G4cout << " " << invalidIdentifier << G4endl;
207 }
208 } else if (id == invalidIdentifier) {
209 G4cout
210 << this->GetName()
211 << " WARNING ProcessHits: Something wrong with LArG4Identifier: "
212 << (std::string)larhit.id << G4endl << " " << id
213 << G4endl << " " << id.getString() << G4endl
214 << " "
215 << a_step->GetPreStepPoint()
216 ->GetPhysicalVolume()
217 ->GetLogicalVolume()
218 ->GetName()
219 << G4endl
220 << " numberOfProcessedHits: " << numberOfProcessedHits
221 << G4endl;
222 G4cout << " " << invalidIdentifier << G4endl;
223 }
224 // need to get the cell information
225 if (numberOfProcessedHits > 1) {
226 if (!m_larEmID->is_em_barrel(id)) {
227 // It didn't seem to happen outside em_barrel, so flag up if it
228 // does:
229 G4cout << this->GetName()
230 << " WARNING ProcessHits: Outside LAr barrel, but "
231 "numberOfProcessedHits="
232 << numberOfProcessedHits
233 << ", LArG4Identifier: " << (std::string)larhit.id << G4endl;
234 } else {
235 if (m_config.shift_lar_subhit) {
236 // find subhit with largest energy
237 if (m_config.verboseLevel > 9) {
238 G4cout << this->GetName()
239 << " VERBOSE ProcessHits: shifting subhits: largest "
240 "energy subhit index is "
241 << maxSubHitEnergyindex << " E: " << maxSubHitEnergy
242 << " identifier: " << maxEnergyIdentifier.getString()
243 << G4endl;
244 }
245 // from identifier
246 const CaloDetDescrElement* thiscell =
247 m_calo_dd_man.get()->get_element(id);
248 if (!maxEnergyCell) {
249 // How often does this happen? Do not shift.
250 G4cout << this->GetName()
251 << " WARNING ProcessHits: maxEnergyCell failed: "
252 << maxEnergyIdentifier.getString() << G4endl
253 << " "
254 << m_calo_dd_man.get()->get_element(id)->getSampling()
255 << G4endl << " " << originalStepPosition.eta()
256 << G4endl << " " << originalStepPosition.phi()
257 << G4endl;
258 stepPosition = originalStepPosition;
259 } else if (maxEnergyCell == thiscell) {
260 // The cells match, so do not shift this hit.
261 if (m_config.verboseLevel > 9) {
262 G4cout << this->GetName()
263 << " VERBOSE ProcessHits: Original step position: "
264 << originalStepPosition.x() << " "
265 << originalStepPosition.y() << " "
266 << originalStepPosition.z() << G4endl << " "
267 << "This cell: " << thiscell->x() << " "
268 << thiscell->y() << " " << thiscell->z() << G4endl
269 << " " << "No shift" << G4endl;
270 }
271 stepPosition = originalStepPosition;
272 } else {
273 // the two cells do not match => shift
274 G4ThreeVector diff(thiscell->x() - maxEnergyCell->x(),
275 thiscell->y() - maxEnergyCell->y(),
276 thiscell->z() - maxEnergyCell->z());
277 stepPosition = originalStepPosition + diff;
278 if (m_config.verboseLevel > 9) {
279 const CaloDetDescrElement* bestcell =
280 m_calo_dd_man.get()->get_element(
281 m_calo_dd_man.get()->get_element(id)->getSampling(),
282 originalStepPosition.eta(),
283 originalStepPosition.phi());
284 if (bestcell) [[likely]] {
285 G4cout << this->GetName()
286 << " VERBOSE ProcessHits: Original step position: "
287 << originalStepPosition.x() << " "
288 << originalStepPosition.y() << " "
289 << originalStepPosition.z() << G4endl << " "
290 << "This cell: " << thiscell->x() << " "
291 << thiscell->y() << " " << thiscell->z() << G4endl
292 << " "
293 << "Highest E cell: " << maxEnergyCell->x() << " "
294 << maxEnergyCell->y() << " " << maxEnergyCell->z()
295 << G4endl << " "
296 << "(Best cell: " << bestcell->x() << " "
297 << bestcell->y() << " " << bestcell->z() << ")"
298 << G4endl << " "
299 << "Shifted step position: " << stepPosition.x() << " "
300 << stepPosition.y() << " " << stepPosition.z()
301 << G4endl;
302 }
303 }
304 }
305 }
306 }
307 }
308 // Finalize time for LAr hits?: NO
309 // double time = larhit.energy)==0 ? 0. : (double)
310 // larhit.time/larhit.energy/CLHEP::ns;
311 double time = larhit.time;
312 double energy = larhit.energy / CLHEP::MeV;
313 // Drop any hits that don't have a good identifier attached
314 if (!m_calo_dd_man.get()->get_element(id)) {
315 if (m_config.verboseLevel > 4) {
316 G4cout << this->GetName()
317 << " DEBUG update_map: bad identifier: " << id.getString()
318 << " skipping this hit." << G4endl;
319 }
320 continue;
321 }
322 this->update_map(
323 stepPosition, id, energy, time, true,
324 numberOfProcessedHits); // store numberOfProcessedHits as info
325 } // numberOfProcessedHits
326 } // istep
327 if (madeSubSteps) {
328 // only delete steps when doing substeps. Do not delete the original
329 // G4Step!
330 while (!steps.empty()) {
331 delete steps.back();
332 steps.pop_back();
333 }
334 }
335 if (m_config.verboseLevel > 4) {
336 G4cout << "LArFCS_StepInfoSD::ProcessHits(): Etotal substeps=" << et_all
337 << G4endl << G4endl << G4endl;
338 }
339 }
340
341 return result;
342}
343
345 Identifier id;
346 if (a_ident[0] == 4) {
347 // is LAr
348 if (a_ident[1] == 1) {
349 // is LAr EM
350 try {
351 id = m_larEmID->channel_id(a_ident[2], // barrel_ec
352 a_ident[3], // sampling
353 a_ident[4], // region
354 a_ident[5], // eta
355 a_ident[6]); // phi
356 } catch (LArID_Exception& e) {
357 G4ExceptionDescription description;
358 description << "ConvertID: LArEM_ID error code " << e.code() << " "
359 << (std::string)e;
360 G4Exception("LArFCS_StepInfoSD", "ConvertIDEM", FatalException,
362 abort();
363 }
364 } else if (a_ident[1] == 2) {
365 // is EM HEC
366 try {
367 id = m_larHecID->channel_id(a_ident[2], // zSide
368 a_ident[3], // sampling
369 a_ident[4], // region
370 a_ident[5], // eta
371 a_ident[6]); // phi
372 } catch (LArID_Exception& e) {
373 G4ExceptionDescription description;
374 description << "ConvertID: LArHEC_ID error code " << e.code() << " "
375 << (std::string)e;
376 G4Exception("LArFCS_StepInfoSD", "ConvertIDHEC", FatalException,
378 abort();
379 }
380 } else if (a_ident[1] == 3) {
381 // FCAL
382 if (a_ident[3] > 0) {
383 // is EM FCAL
384 try {
385 id = m_larFcalID->channel_id(a_ident[2], // zSide
386 a_ident[3], // sampling
387 a_ident[4], // eta
388 a_ident[5]); // phi
389 } catch (LArID_Exception& e) {
390 G4ExceptionDescription description;
391 description << "ConvertID: LArFCAL_ID error code " << e.code() << " "
392 << (std::string)e;
393 G4Exception("LArFCS_StepInfoSD", "ConvertIDFCAL", FatalException,
395 abort();
396 }
397 } else {
398 G4ExceptionDescription description;
399 description << "ConvertID: Unsupported ID. ";
400 G4Exception("LArFCS_StepInfoSD", "ConvertIDMiniFCAL", FatalException,
402 abort();
403 }
404 }
405 }
406 return id;
407}
Scalar mag() const
mag method
Definition of CaloDetDescrManager.
float et(const xAOD::jFexSRJetRoI *j)
void diff(const Jet &rJet1, const Jet &rJet2, std::map< std::string, double > varDiff)
Difference between jets - Non-Class function required by trigger.
Definition Jet.cxx:631
This class groups all DetDescr information related to a CaloCell.
FCS_Param::Config m_config
void update_map(const CLHEP::Hep3Vector &l_vec, const Identifier &l_identifier, double l_energy, double l_time, bool l_valid, int l_detector)
CxxUtils::CachedPointer< const CaloDetDescrManager > m_calo_dd_man
const LArHEC_ID * m_larHecID
FCS_StepInfoSD(G4String a_name, const FCS_Param::Config &config)
Constructor.
const LArFCAL_ID * m_larFcalID
const LArEM_ID * m_larEmID
Pointers to the identifier helpers.
void getCaloDDManager()
Keep a map instead of trying to keep the full vector.
std::string getString() const
Provide a string form of the identifier - hexadecimal.
LArFCS_StepInfoSD(G4String a_name, const FCS_Param::Config &config)
Constructor.
G4bool ProcessHits(G4Step *a_step, G4TouchableHistory *) override
Main processing method.
virtual ~LArFCS_StepInfoSD()
Destructor.
ILArCalculatorSvc * m_calculator
Member variable - the calculator we'll use.
Identifier ConvertID(const LArG4Identifier &a_ident) const
Helper function for making "real" identifiers from LArG4Identifiers.
Exception class for LAr Identifiers.
std::string description
glabal timer - how long have I taken so far?
Definition hcg.cxx:91
STL namespace.
#define likely(x)
Extra patterns decribing particle interation process.