ATLAS Offline Software
FastCaloSim.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 // Header include
6 #include "FastCaloSim.h"
7 
8 // FastCaloSim includes
12 
13 // Geant4 particle includes
14 #include "G4Gamma.hh"
15 #include "G4Electron.hh"
16 #include "G4Positron.hh"
17 #include "G4PionPlus.hh"
18 #include "G4PionMinus.hh"
19 
20 //Geant4
21 #include "G4ParticleTable.hh"
22 #include "Randomize.hh"
23 
24 // HepMCHelpers include
26 
27 // G4 sensitive detector includes
28 #include "G4SDManager.hh"
30 
31 #undef FCS_DEBUG
32 
33 FastCaloSim::FastCaloSim(const std::string& name,
34  G4Region* region,
35  const PublicToolHandle<IFastCaloSimCaloTransportation>& FastCaloSimCaloTransportation,
36  const PublicToolHandle<IFastCaloSimCaloExtrapolation>& FastCaloSimCaloExtrapolation,
37  const PublicToolHandle<IG4CaloTransportTool>& G4CaloTransportTool,
38  const PublicToolHandle<IPunchThroughSimWrapper>& PunchThroughSimWrapper,
39  const ServiceHandle<ISF::IFastCaloSimParamSvc>& FastCaloSimSvc,
40  const std::string& CaloCellContainerSDName,
41  bool doG4Transport,
42  bool doPunchThrough,
44 
45 : G4VFastSimulationModel(name, region),
46  m_FastCaloSimCaloTransportation(FastCaloSimCaloTransportation),
47  m_FastCaloSimCaloExtrapolation(FastCaloSimCaloExtrapolation),
48  m_G4CaloTransportTool(G4CaloTransportTool),
49  m_PunchThroughSimWrapper(PunchThroughSimWrapper),
50  m_FastCaloSimSvc(FastCaloSimSvc),
51  m_CaloCellContainerSDName(CaloCellContainerSDName),
52  m_doG4Transport(doG4Transport),
53  m_doPunchThrough(doPunchThrough),
54  m_FastCaloSimTool(FastCaloSimTool)
55 {
56 }
57 
58 
59 G4bool FastCaloSim::IsApplicable(const G4ParticleDefinition& particleType)
60 {
61  // Check whether we can simulate the particle with FastCaloSim
62  bool isPhoton = &particleType == G4Gamma::GammaDefinition();
63  bool isElectron = &particleType == G4Electron::ElectronDefinition();
64  bool isPositron = &particleType == G4Positron::PositronDefinition();
65  bool isHadron = MC::isHadron(particleType.GetPDGEncoding());
66 
67  // FastCaloSim is applicable if it is photon, electron, positron or any hadron
68  bool isApplicable = isPhoton || isElectron || isPositron || isHadron;
69 
70  #ifdef FCS_DEBUG
71  const std::string pName = particleType.GetParticleName();
72  G4cout<< "[FastCaloSim::IsApplicable] Got " << pName <<G4endl;
73  if(isApplicable) G4cout<<"[FastCaloSim::IsApplicable] APPLICABLE"<<G4endl;
74  else G4cout<<"[FastCaloSim::IsApplicable] NOT APPLICABLE"<<G4endl;
75  #endif
76 
77 
78  return isApplicable;
79 }
80 
81 G4bool FastCaloSim::ModelTrigger(const G4FastTrack& fastTrack)
82 {
83 
84  #ifdef FCS_DEBUG
85  G4cout<<"[FastCaloSim::ModelTrigger] Got particle with " <<"\n"
86  <<" pdg=" <<fastTrack.GetPrimaryTrack() -> GetDefinition()->GetPDGEncoding() <<"\n"
87  <<" Ekin="<<fastTrack.GetPrimaryTrack() -> GetKineticEnergy() <<"\n"
88  <<" p=" <<fastTrack.GetPrimaryTrack() -> GetMomentum().mag() <<"\n"
89  <<" x=" <<fastTrack.GetPrimaryTrack() -> GetPosition().x() <<"\n"
90  <<" y=" <<fastTrack.GetPrimaryTrack() -> GetPosition().y() <<"\n"
91  <<" z=" <<fastTrack.GetPrimaryTrack() -> GetPosition().z() <<"\n"
92  <<" r=" <<fastTrack.GetPrimaryTrack() -> GetPosition().perp() <<"\n"
93  <<" eta=" <<fastTrack.GetPrimaryTrack() -> GetMomentum().eta() <<"\n"
94  <<" phi=" <<fastTrack.GetPrimaryTrack() -> GetMomentum().phi() <<"\n"
95  <<G4endl;
96  #endif
97 
98 
99  // Simulate particles below 50 keV with Geant4 to have same config as ISF implementation
100  if (fastTrack.GetPrimaryTrack() -> GetKineticEnergy() < 0.05) {
101  #ifdef FCS_DEBUG
102  G4cout<<"[FastCaloSim::ModelTrigger] Particle below 50 keV threshold. Passing to G4. "<<G4endl;
103  #endif
104  return false;
105  }
106 
107  // Check if triggered particle is really on the ID-Calo (parametrization) boundary
108  if (!passedIDCaloBoundary(fastTrack)) {
109  #ifdef FCS_DEBUG
110  G4cout<<"[FastCaloSim::ModelTrigger] Particle failed passedIDCaloBoundary z="<<fastTrack.GetPrimaryTrack() -> GetPosition().z()<<" r="<<fastTrack.GetPrimaryTrack() -> GetPosition().perp()<<G4endl;
111  #endif
112  return false;
113  }
114 
115  // Set minimum kinetic energy of pions and other hadrons required to be passed to FastCaloSim
116  float minEkinPions = 200;
117  float minEkinOtherHadrons = 400;
118 
119  // Get particle definition
120  const G4ParticleDefinition * G4Particle = fastTrack.GetPrimaryTrack() -> GetDefinition();
121  // Get particle kinetic energy
122  const float Ekin = fastTrack.GetPrimaryTrack() -> GetKineticEnergy();
123 
124  // Check particle type
125  bool isPhoton = G4Particle == G4Gamma::Definition();
126  bool isElectron = G4Particle == G4Electron::Definition();
127  bool isPositron = G4Particle == G4Positron::Definition();
128  bool isPionPlus = G4Particle == G4PionPlus::Definition();
129  bool isPionMinus = G4Particle == G4PionMinus::Definition();
130 
131  // Pass all photons, electrons and positrons to FastCaloSim
132  if (isPhoton || isElectron || isPositron){
133  #ifdef FCS_DEBUG
134  G4cout<<"[FastCaloSim::ModelTrigger] Model triggered"<<G4endl;
135  #endif
136  return true;
137  }
138 
139  // Require minimum kinetic energy of pions needed to be passed to FastCaloSim
140  if (isPionPlus || isPionMinus){
141  bool passMinEkinPions = Ekin > minEkinPions;
142 
143  #ifdef FCS_DEBUG
144  if(passMinEkinPions) G4cout<<"[FastCaloSim::ModelTrigger] Model triggered"<<G4endl;
145  else G4cout<<"[FastCaloSim::ModelTrigger] Pion with Ekin="<<Ekin<<" below the minimum "<<minEkinPions<<" MeV threshold. Model not triggered."<<G4endl;
146  #endif
147 
148  return passMinEkinPions;
149 
150  }
151 
152  // Require minimum kinetic energy of other hadrons needed to be passed to FastCaloSim
153  bool passMinEkinOtherHadrons = Ekin > minEkinOtherHadrons;
154 
155  #ifdef FCS_DEBUG
156  if(passMinEkinOtherHadrons) G4cout<<"[FastCaloSim::ModelTrigger] Model triggered"<<G4endl;
157  else G4cout<<"[FastCaloSim::ModelTrigger] Other hadron with Ekin="<<Ekin<<" below the minimum "<<minEkinOtherHadrons<<" MeV threshold. Model not triggered."<<G4endl;
158  #endif
159 
160  return passMinEkinOtherHadrons;
161 }
162 
163 void FastCaloSim::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep)
164 {
165 
166  TFCSSimulationState simState(G4Random::getTheEngine());
167  TFCSTruthState truthState;
168  TFCSExtrapolationState extrapolState;
169 
170  // Get Geant4 primary track
171  const G4Track * G4PrimaryTrack = fastTrack.GetPrimaryTrack();
172  // Get Geant4 particle definition
173  const G4ParticleDefinition * G4Particle = G4PrimaryTrack -> GetDefinition();
174  // Get Geant4 particle pdgID
175  signed int pdgID = G4Particle -> GetPDGEncoding();
176 
177  // Do not simulate particles below 10 MeV
178  if(G4PrimaryTrack -> GetKineticEnergy() < 10){
179  #ifdef FCS_DEBUG
180  G4cout<<"[FastCaloSim::DoIt] Skipping particle with Ekin: " << G4PrimaryTrack -> GetKineticEnergy() <<" MeV. Below the 10 MeV threshold"<<G4endl;
181  #endif
182  fastStep.KillPrimaryTrack();
183  return;
184  }
185  // Decide on which FastCaloSim parametrization to use (electron, photon or pion)
186  if(G4Particle == G4Electron::Definition() || G4Particle == G4Positron::Definition() || G4Particle == G4Gamma::Definition())
187  {
188  // Use egamma parametrization for simulation of electrons and photons
189  truthState.set_pdgid(pdgID);
190  }
191  else{
192  // Use pion parametrization for simulation of all hadrons
193  truthState.set_pdgid(G4PionPlus::Definition() -> GetPDGEncoding());
194  }
195 
196  // Set the kinematics of the FastCaloSim truth state
197  truthState.SetPtEtaPhiM(G4PrimaryTrack -> GetMomentum().perp(),
198  G4PrimaryTrack -> GetMomentum().eta(),
199  G4PrimaryTrack -> GetMomentum().phi(),
200  G4Particle -> GetPDGMass());
201 
202  // Set the vertex of the FastCaloSim truth state
203  truthState.set_vertex(G4PrimaryTrack -> GetPosition().x(),
204  G4PrimaryTrack -> GetPosition().y(),
205  G4PrimaryTrack -> GetPosition().z());
206 
207 
208  /* For anti protons and anti-neutrons the kinetic energy should be
209  calculated as Ekin = E() + M() instead of E() - M() this is
210  achieved by setting an Ekin offset of 2*M() to the truth state */
211  if(pdgID == -2212 || pdgID == -2112) truthState.set_Ekin_off(2 * G4Particle -> GetPDGMass());
212 
213 
214  // Perform particle transportation through calorimeter system eiher with ATLAS tracking tools or with Geant4
215  std::vector<G4FieldTrack> caloSteps = m_doG4Transport ? m_G4CaloTransportTool -> transport(*G4PrimaryTrack)
216  : m_FastCaloSimCaloTransportation -> transport(&truthState, false);
217 
218  // Extrapolate transported stepos to ID-Calo boundary and all layers of the calorimeter system
219  m_FastCaloSimCaloExtrapolation->extrapolate(extrapolState, &truthState, caloSteps);
220 
221  // Do not simulate further if extrapolation to ID - Calo boundary fails
222  if(extrapolState.IDCaloBoundary_eta() == -999){
223  #ifdef FCS_DEBUG
224  G4cout<<"[FastCaloSim::DoIt] Killing particle as extrapolation failed"<<G4endl;
225  #endif
226  fastStep.KillPrimaryTrack();
227  return;
228  }
229  // Perform the actual FastCaloSim simulation
230  if(m_FastCaloSimSvc->simulate(simState, &truthState, &extrapolState).isFailure()){
231  G4Exception("FastCaloSimSvc", "FailedSimulationCall", FatalException, "FastCaloSimSvc: Simulation call failed.");
232  abort();
233  }
234 
235  #ifdef FCS_DEBUG
236  G4cout<<"[FastCaloSim::DoIt] pdgID of G4PrimaryTrack: " << pdgID << G4endl;
237  G4cout<<"[FastCaloSim::DoIt] Energy returned: " << simState.E() << G4endl;
238  G4cout<<"[FastCaloSim::DoIt] Energy fraction for layer: " << G4endl;
239  for (int s = 0; s < 24; s++) G4cout<<"[FastCaloSim::DoIt] Sampling " << s << " energy " << simState.E(s) << G4endl;
240  #endif
241 
242  // Retrieve the associated CaloCellContainer sensitive detector
243  CaloCellContainerSD * caloCellContainerSD = getCaloCellContainerSD();
244  // Record the cells
245  caloCellContainerSD->recordCells(simState);
246 
247  // Do punchthrough here (secondaries), after the main simulation
248  if (m_doPunchThrough){
249  // necessary for determining particle type and properties (mass etc)
250  G4ParticleTable *ptable = G4ParticleTable::GetParticleTable();
251 
252  // Get simulated energy
253  const double simE = simState.E();
254 
255  // Get energy fraction in layers into vector
256  std::vector<double> simEfrac;
257  for (unsigned int i = 0; i < 24; i++){simEfrac.push_back(simState.Efrac(i));}
258 
259  // run actual method (no return, it will do fastStep.CreateSecondaryTrack(...) under the hood)
260  m_PunchThroughSimWrapper->DoPunchThroughSim(*ptable, G4Random::getTheEngine(), simE, simEfrac, fastTrack, fastStep);
261  }
262 
263  // Clean up the auxiliary info from the simulation state
264  simState.DoAuxInfoCleanup();
265 
266  // Finally kill the primary track after all simulation steps done
267  fastStep.KillPrimaryTrack();
268  fastStep.ProposePrimaryTrackPathLength(0.0);
269 }
270 
272 
273  G4SDManager *sdm = G4SDManager::GetSDMpointer();
274  G4VSensitiveDetector * vsd = sdm->FindSensitiveDetector((G4String)m_CaloCellContainerSDName);
275 
276  if (!vsd){
277  G4Exception("FastCaloSimSvc", "FailedFindSensitiveDetector", FatalException, "FastCaloSimSvc: Failed getting CaloCellContainer SD.");
278  abort();
279  }
280  // Cast G4VSensitiveDetector to CaloCellContainerSD
281  CaloCellContainerSD * caloCellContainerSD = dynamic_cast<CaloCellContainerSD*>(vsd);
282 
283  if (!caloCellContainerSD){
284  G4Exception("FastCaloSimSvc", "FailedCastSensitiveDetector", FatalException, "FastCaloSimSvc: Failed casting G4VSensitiveDetector.");
285  abort();
286  }
287  return caloCellContainerSD;
288 }
289 
290 
291 G4bool FastCaloSim::passedIDCaloBoundary(const G4FastTrack& fastTrack){
292 
293 
294  /* Method checks if particle has crossed the ID-Calo boundary, defined using three cylinders with pairs of r/z values
295  We also perform a directional check to make sure that the particle does not originate from any backscatter from the MS or CALO */
296 
297  // NOTE:
298  // For the inner beam pipe section, innerBeamPipeZ = 4185 corresponds to original AFII boundary, while
299  // innerBeamPipeZ = 4587 corresponds to the CALO::CALO boundary which should be used instead as there
300  // is no triggering volume at this point and we will trigger slightly later than the AFII boundary, so
301  // passedIDCaloBoundary would fail and we would simulate this part with Geant4
302 
303  // Barrel
304  const float barrelR = 1148;
305  const float barrelZ = 3549.5;
306  // Inner beam pipe section
307  const float innerBeamPipeR = 120;
308  const float innerBeamPipeZ = 4587;
309  // Outer beam pipe section
310  const float outerBeamPipeR = 41;
311  const float outerBeamPipeZ = 6783;
312 
313  // Get particle position
314  const G4ThreeVector particlePosition = fastTrack.GetPrimaryTrack() -> GetPosition();
315  // Get r and z values of position
316  const float r = particlePosition.perp();
317  const float z = particlePosition.z();
318  // Get direction of particle
319  const G4ThreeVector particleDirection = fastTrack.GetPrimaryTrack() -> GetMomentum();
320  // Construct the helper line to decide if particle comes from ID or is backscatter from CALO
321  const G4ThreeVector helperLine = particlePosition + particleDirection;
322 
323  // Set 5cm trigger tolerance, i.e. the 'width' of the trigger boundary
324  // be careful not to set this too low, else Geant4 might overjump your trigger boundary and simulate everything with G4
325  const float triggerTolerance = 50;
326 
327  // Barrel boundary (horizontal surfaces)
328  if (std::abs(z) <= barrelZ + triggerTolerance) {
329  if (r >= barrelR && r < barrelR + triggerTolerance) return helperLine.perp() >= barrelR;
330  }
331  // Beam pipe boundary (horizontal surfaces)
332  if (std::abs(z) >= barrelZ && std::abs(z) <= innerBeamPipeZ + triggerTolerance){
333  if (r >= innerBeamPipeR && r < innerBeamPipeR + triggerTolerance) return helperLine.perp() >= innerBeamPipeR;
334  }
335  if (std::abs(z) >= innerBeamPipeZ && std::abs(z) <= outerBeamPipeZ){
336  if (r >= outerBeamPipeR && r < outerBeamPipeR + triggerTolerance) return helperLine.perp() >= outerBeamPipeR;
337  }
338  // Barrel boundary (vertical surfaces)
339  if (r >= outerBeamPipeR && r<= innerBeamPipeR){
340  if (std::abs(z) >= innerBeamPipeZ && std::abs(z) < innerBeamPipeZ + triggerTolerance) return std::abs(helperLine.z()) >= innerBeamPipeZ;
341  }
342  // Beam pipe boundary (vertical surfaces)
343  if (r >= innerBeamPipeR && r <= barrelR){
344  if (std::abs(z) >= barrelZ && std::abs(z) < barrelZ + triggerTolerance) return std::abs(helperLine.z()) >= barrelZ;
345  }
346 
347  return false;
348 
349 }
350 
FastCaloSimCaloExtrapolation
Definition: FastCaloSimCaloExtrapolation.h:32
CaloCellContainerSD.h
beamspotman.r
def r
Definition: beamspotman.py:672
TFCSSimulationState::DoAuxInfoCleanup
void DoAuxInfoCleanup()
Definition: TFCSSimulationState.cxx:79
phi
Scalar phi() const
phi method
Definition: AmgMatrixBasePlugin.h:67
perp
Scalar perp() const
perp method - perpenticular length
Definition: AmgMatrixBasePlugin.h:44
TFCSTruthState::set_Ekin_off
void set_Ekin_off(double val)
Definition: TFCSTruthState.h:23
eta
Scalar eta() const
pseudorapidity method
Definition: AmgMatrixBasePlugin.h:83
CaloCellContainerSD::recordCells
void recordCells(TFCSSimulationState &)
Definition: CaloCellContainerSD.cxx:65
TFCSSimulationState::Efrac
double Efrac(int sample) const
Definition: TFCSSimulationState.h:44
FastCaloSim::m_FastCaloSimSvc
ServiceHandle< ISF::IFastCaloSimParamSvc > m_FastCaloSimSvc
Definition: FastCaloSim.h:73
TFCSExtrapolationState
Definition: TFCSExtrapolationState.h:13
x
#define x
particleType
Definition: particleType.h:29
FastCaloSim::FastCaloSim
FastCaloSim(const std::string &name, G4Region *region, const PublicToolHandle< IFastCaloSimCaloTransportation > &FastCaloSimCaloTransportation, const PublicToolHandle< IFastCaloSimCaloExtrapolation > &FastCaloSimCaloExtrapolation, const PublicToolHandle< IG4CaloTransportTool > &G4CaloTransportTool, const PublicToolHandle< IPunchThroughSimWrapper > &PunchThroughSimWrapper, const ServiceHandle< ISF::IFastCaloSimParamSvc > &FastCaloSimSvc, const std::string &CaloCellContainerSDName, bool doG4Transport, bool doPunchThrough, FastCaloSimTool *FastCaloSimTool)
Definition: FastCaloSim.cxx:33
CaloCellContainerSD
Definition: CaloCellContainerSD.h:22
FastCaloSim::getCaloCellContainerSD
CaloCellContainerSD * getCaloCellContainerSD()
Retrieves the associated sensitive detector responsible for writing out the CaloCellContainer.
Definition: FastCaloSim.cxx:271
FastCaloSim::m_FastCaloSimCaloTransportation
PublicToolHandle< IFastCaloSimCaloTransportation > m_FastCaloSimCaloTransportation
Definition: FastCaloSim.h:64
FastCaloSim.h
FastCaloSim::m_doPunchThrough
bool m_doPunchThrough
Definition: FastCaloSim.h:80
FastCaloSim::m_CaloCellContainerSDName
std::string m_CaloCellContainerSDName
Definition: FastCaloSim.h:75
lumiFormat.i
int i
Definition: lumiFormat.py:85
z
#define z
xAOD::EgammaHelpers::isElectron
bool isElectron(const xAOD::Egamma *eg)
is the object an electron (not Fwd)
Definition: EgammaxAODHelpers.cxx:12
TFCSExtrapolationState::IDCaloBoundary_eta
double IDCaloBoundary_eta() const
Definition: TFCSExtrapolationState.h:63
FastCaloSim::m_PunchThroughSimWrapper
PublicToolHandle< IPunchThroughSimWrapper > m_PunchThroughSimWrapper
Definition: FastCaloSim.h:70
FastCaloSim::m_FastCaloSimCaloExtrapolation
PublicToolHandle< IFastCaloSimCaloExtrapolation > m_FastCaloSimCaloExtrapolation
Definition: FastCaloSim.h:66
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:240
isHadron
bool isHadron(const T &p)
Definition: AtlasPID.h:351
PunchThroughSimWrapper
Class to wrap PunchThrough simulation inside FastCaloSim; Runs both PunchThroughG4Classifier and Punc...
Definition: PunchThroughSimWrapper.h:32
FastCaloSim::ModelTrigger
G4bool ModelTrigger(const G4FastTrack &) override final
Determines the applicability of the fast sim model to this particular track.
Definition: FastCaloSim.cxx:81
FastCaloSim::m_doG4Transport
bool m_doG4Transport
Definition: FastCaloSim.h:77
TFCSTruthState.h
y
#define y
TFCSExtrapolationState.h
xAOD::EgammaHelpers::isPhoton
bool isPhoton(const xAOD::Egamma *eg)
is the object a photon
Definition: EgammaxAODHelpers.cxx:21
FastCaloSimCaloTransportation
Definition: FastCaloSimCaloTransportation.h:21
TFCSTruthState::set_vertex
void set_vertex(const TLorentzVector &val)
Definition: TFCSTruthState.h:19
python.SystemOfUnits.s
float s
Definition: SystemOfUnits.py:147
FastCaloSim::DoIt
void DoIt(const G4FastTrack &, G4FastStep &) override final
Definition: FastCaloSim.cxx:163
TFCSSimulationState::E
double E() const
Definition: TFCSSimulationState.h:42
TFCSSimulationState.h
FastCaloSim::m_G4CaloTransportTool
PublicToolHandle< IG4CaloTransportTool > m_G4CaloTransportTool
Definition: FastCaloSim.h:68
FastCaloSim::passedIDCaloBoundary
G4bool passedIDCaloBoundary(const G4FastTrack &fastTrack)
Check if the particle is located at the proper ID-Calo parametrization boundary and is travelling out...
Definition: FastCaloSim.cxx:291
G4CaloTransportTool
A tool which transports particles through the Geant4 geometry.
Definition: G4CaloTransportTool.h:25
TFCSTruthState
Definition: TFCSTruthState.h:13
TFCSSimulationState
Definition: TFCSSimulationState.h:32
FastCaloSim::IsApplicable
G4bool IsApplicable(const G4ParticleDefinition &) override final
Definition: FastCaloSim.cxx:59
FastCaloSimTool
Definition: FastCaloSimTool.h:24
HepMCHelpers.h
TFCSTruthState::set_pdgid
void set_pdgid(int val)
Definition: TFCSTruthState.h:18
ServiceHandle< ISF::IFastCaloSimParamSvc >