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