ATLAS Offline Software
Loading...
Searching...
No Matches
FastCaloSim Class Reference

#include <FastCaloSim.h>

Inheritance diagram for FastCaloSim:
Collaboration diagram for FastCaloSim:

Public Member Functions

 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)
 ~FastCaloSim ()
G4bool IsApplicable (const G4ParticleDefinition &) override final
void DoIt (const G4FastTrack &, G4FastStep &) override final
G4bool ModelTrigger (const G4FastTrack &) override final
 Determines the applicability of the fast sim model to this particular track.
CaloCellContainerSDgetCaloCellContainerSD ()
 Retrieves the associated sensitive detector responsible for writing out the CaloCellContainer.
G4bool passedIDCaloBoundary (const G4FastTrack &fastTrack)
 Check if the particle is located at the proper ID-Calo parametrization boundary and is travelling outwards from the ID to the CALO.

Private Attributes

PublicToolHandle< IFastCaloSimCaloTransportationm_FastCaloSimCaloTransportation
PublicToolHandle< IFastCaloSimCaloExtrapolationm_FastCaloSimCaloExtrapolation
PublicToolHandle< IG4CaloTransportToolm_G4CaloTransportTool
PublicToolHandle< IPunchThroughSimWrapperm_PunchThroughSimWrapper
ServiceHandle< ISF::IFastCaloSimParamSvcm_FastCaloSimSvc
std::string m_CaloCellContainerSDName
bool m_doG4Transport
bool m_doPunchThrough

Detailed Description

Definition at line 31 of file FastCaloSim.h.

Constructor & Destructor Documentation

◆ 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 at line 36 of file FastCaloSim.cxx.

48: G4VFastSimulationModel(name, region),
49 m_FastCaloSimCaloTransportation(FastCaloSimCaloTransportation),
50 m_FastCaloSimCaloExtrapolation(FastCaloSimCaloExtrapolation),
51 m_G4CaloTransportTool(G4CaloTransportTool),
52 m_PunchThroughSimWrapper(PunchThroughSimWrapper),
53 m_FastCaloSimSvc(FastCaloSimSvc),
54 m_CaloCellContainerSDName(CaloCellContainerSDName),
55 m_doG4Transport(doG4Transport),
56 m_doPunchThrough(doPunchThrough)
57{
58}
std::string m_CaloCellContainerSDName
Definition FastCaloSim.h:75
PublicToolHandle< IFastCaloSimCaloExtrapolation > m_FastCaloSimCaloExtrapolation
Definition FastCaloSim.h:66
PublicToolHandle< IG4CaloTransportTool > m_G4CaloTransportTool
Definition FastCaloSim.h:68
bool m_doG4Transport
Definition FastCaloSim.h:77
ServiceHandle< ISF::IFastCaloSimParamSvc > m_FastCaloSimSvc
Definition FastCaloSim.h:73
bool m_doPunchThrough
Definition FastCaloSim.h:80
PublicToolHandle< IFastCaloSimCaloTransportation > m_FastCaloSimCaloTransportation
Definition FastCaloSim.h:64
PublicToolHandle< IPunchThroughSimWrapper > m_PunchThroughSimWrapper
Definition FastCaloSim.h:70

◆ ~FastCaloSim()

FastCaloSim::~FastCaloSim ( )
inline

Definition at line 46 of file FastCaloSim.h.

46{}

Member Function Documentation

◆ DoIt()

void FastCaloSim::DoIt ( const G4FastTrack & fastTrack,
G4FastStep & fastStep )
finaloverride

Definition at line 165 of file FastCaloSim.cxx.

166{
167 if (auto* eventInfo = AtlasG4EventUserInfo::GetEventUserInfo()) {
168 Gaudi::Hive::setCurrentContext(eventInfo->GetEventContext());
169 }
170
171 TFCSSimulationState simState(G4Random::getTheEngine());
172 TFCSTruthState truthState;
173 TFCSExtrapolationState extrapolState;
174
175 // Get Geant4 primary track
176 const G4Track * G4PrimaryTrack = fastTrack.GetPrimaryTrack();
177 // Get Geant4 particle definition
178 const G4ParticleDefinition * G4Particle = G4PrimaryTrack -> GetDefinition();
179 // Get Geant4 particle pdgID
180 signed int pdgID = G4Particle -> GetPDGEncoding();
181
182 // Do not simulate particles below 10 MeV
183 if(G4PrimaryTrack -> GetKineticEnergy() < 10){
184 #ifdef FCS_DEBUG
185 G4cout<<"[FastCaloSim::DoIt] Skipping particle with Ekin: " << G4PrimaryTrack -> GetKineticEnergy() <<" MeV. Below the 10 MeV threshold"<<G4endl;
186 #endif
187 fastStep.KillPrimaryTrack();
188 return;
189 }
190 // Decide on which FastCaloSim parametrization to use (electron, photon or pion)
191 if(G4Particle == G4Electron::Definition() || G4Particle == G4Positron::Definition() || G4Particle == G4Gamma::Definition())
192 {
193 // Use egamma parametrization for simulation of electrons and photons
194 truthState.set_pdgid(pdgID);
195 }
196 else{
197 // Use pion parametrization for simulation of all hadrons
198 truthState.set_pdgid(G4PionPlus::Definition() -> GetPDGEncoding());
199 }
200
201 // Set the kinematics of the FastCaloSim truth state
202 truthState.SetPtEtaPhiM(G4PrimaryTrack -> GetMomentum().perp(),
203 G4PrimaryTrack -> GetMomentum().eta(),
204 G4PrimaryTrack -> GetMomentum().phi(),
205 G4Particle -> GetPDGMass());
206
207 // Set the vertex of the FastCaloSim truth state
208 truthState.set_vertex(G4PrimaryTrack -> GetPosition().x(),
209 G4PrimaryTrack -> GetPosition().y(),
210 G4PrimaryTrack -> GetPosition().z());
211
212
213 /* For anti protons and anti-neutrons the kinetic energy should be
214 calculated as Ekin = E() + M() instead of E() - M() this is
215 achieved by setting an Ekin offset of 2*M() to the truth state */
216 if(pdgID == -2212 || pdgID == -2112) truthState.set_Ekin_off(2 * G4Particle -> GetPDGMass());
217
218
219 // Perform particle transportation through calorimeter system eiher with ATLAS tracking tools or with Geant4
220 std::vector<G4FieldTrack> caloSteps = m_doG4Transport ? m_G4CaloTransportTool -> transport(*G4PrimaryTrack)
221 : m_FastCaloSimCaloTransportation -> transport(&truthState, false);
222
223 // Extrapolate transported stepos to ID-Calo boundary and all layers of the calorimeter system
224 m_FastCaloSimCaloExtrapolation->extrapolate(extrapolState, &truthState, caloSteps);
225
226 // Do not simulate further if extrapolation to ID - Calo boundary fails
227 if(extrapolState.IDCaloBoundary_eta() == -999){
228 #ifdef FCS_DEBUG
229 G4cout<<"[FastCaloSim::DoIt] Killing particle as extrapolation failed"<<G4endl;
230 #endif
231 fastStep.KillPrimaryTrack();
232 return;
233 }
234 // Perform the actual FastCaloSim simulation
235 if(m_FastCaloSimSvc->simulate(simState, &truthState, &extrapolState).isFailure()){
236 G4Exception("FastCaloSimSvc", "FailedSimulationCall", FatalException, "FastCaloSimSvc: Simulation call failed.");
237 abort();
238 }
239
240 #ifdef FCS_DEBUG
241 G4cout<<"[FastCaloSim::DoIt] pdgID of G4PrimaryTrack: " << pdgID << G4endl;
242 G4cout<<"[FastCaloSim::DoIt] Energy returned: " << simState.E() << G4endl;
243 G4cout<<"[FastCaloSim::DoIt] Energy fraction for layer: " << G4endl;
244 for (int s = 0; s < 24; s++) G4cout<<"[FastCaloSim::DoIt] Sampling " << s << " energy " << simState.E(s) << G4endl;
245 #endif
246
247 // Retrieve the associated CaloCellContainer sensitive detector
248 CaloCellContainerSD * caloCellContainerSD = getCaloCellContainerSD();
249 // Record the cells
250 caloCellContainerSD->recordCells(simState);
251
252 // Do punchthrough here (secondaries), after the main simulation
253 if (m_doPunchThrough){
254 // necessary for determining particle type and properties (mass etc)
255 G4ParticleTable *ptable = G4ParticleTable::GetParticleTable();
256
257 // Get simulated energy
258 const double simE = simState.E();
259
260 // Get energy fraction in layers into vector
261 std::vector<double> simEfrac;
262 for (unsigned int i = 0; i < 24; i++){simEfrac.push_back(simState.Efrac(i));}
263
264 // run actual method (no return, it will do fastStep.CreateSecondaryTrack(...) under the hood)
265 m_PunchThroughSimWrapper->DoPunchThroughSim(*ptable, G4Random::getTheEngine(), simE, std::move(simEfrac), fastTrack, fastStep);
266 }
267
268 // Clean up the auxiliary info from the simulation state
269 simState.DoAuxInfoCleanup();
270
271 // Finally kill the primary track after all simulation steps done
272 fastStep.KillPrimaryTrack();
273 fastStep.ProposePrimaryTrackPathLength(0.0);
274}
Scalar eta() const
pseudorapidity method
Scalar perp() const
perp method - perpendicular length
Scalar phi() const
phi method
#define y
#define x
#define z
static AtlasG4EventUserInfo * GetEventUserInfo()
void recordCells(TFCSSimulationState &)
CaloCellContainerSD * getCaloCellContainerSD()
Retrieves the associated sensitive detector responsible for writing out the CaloCellContainer.
void set_pdgid(int val)
void set_vertex(const TLorentzVector &val)
void set_Ekin_off(double val)
for(size_t i=0;i< m_blockFillers.size();i++)
Fill one block.

◆ getCaloCellContainerSD()

CaloCellContainerSD * FastCaloSim::getCaloCellContainerSD ( )

Retrieves the associated sensitive detector responsible for writing out the CaloCellContainer.

Definition at line 276 of file FastCaloSim.cxx.

276 {
277
278 G4SDManager *sdm = G4SDManager::GetSDMpointer();
279 G4VSensitiveDetector * vsd = sdm->FindSensitiveDetector((G4String)m_CaloCellContainerSDName);
280
281 if (!vsd){
282 G4Exception("FastCaloSimSvc", "FailedFindSensitiveDetector", FatalException, "FastCaloSimSvc: Failed getting CaloCellContainer SD.");
283 abort();
284 }
285 // Cast G4VSensitiveDetector to CaloCellContainerSD
286 CaloCellContainerSD * caloCellContainerSD = dynamic_cast<CaloCellContainerSD*>(vsd);
287
288 if (!caloCellContainerSD){
289 G4Exception("FastCaloSimSvc", "FailedCastSensitiveDetector", FatalException, "FastCaloSimSvc: Failed casting G4VSensitiveDetector.");
290 abort();
291 }
292 return caloCellContainerSD;
293}

◆ IsApplicable()

G4bool FastCaloSim::IsApplicable ( const G4ParticleDefinition & particleType)
finaloverride

Definition at line 61 of file FastCaloSim.cxx.

62{
63 // Check whether we can simulate the particle with FastCaloSim
64 bool isPhoton = &particleType == G4Gamma::GammaDefinition();
65 bool isElectron = &particleType == G4Electron::ElectronDefinition();
66 bool isPositron = &particleType == G4Positron::PositronDefinition();
67 bool isHadron = MC::isHadron(particleType.GetPDGEncoding());
68
69 // FastCaloSim is applicable if it is photon, electron, positron or any hadron
70 bool isApplicable = isPhoton || isElectron || isPositron || isHadron;
71
72 #ifdef FCS_DEBUG
73 const std::string pName = particleType.GetParticleName();
74 G4cout<< "[FastCaloSim::IsApplicable] Got " << pName <<G4endl;
75 if(isApplicable) G4cout<<"[FastCaloSim::IsApplicable] APPLICABLE"<<G4endl;
76 else G4cout<<"[FastCaloSim::IsApplicable] NOT APPLICABLE"<<G4endl;
77 #endif
78
79
80 return isApplicable;
81}
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
bool isHadron(const T &p)

◆ ModelTrigger()

G4bool FastCaloSim::ModelTrigger ( const G4FastTrack & fastTrack)
finaloverride

Determines the applicability of the fast sim model to this particular track.

Checks that geometric location, energy, and particle type are within bounds. Also checks for containment of the particle's shower within a specific detector region.

Definition at line 83 of file FastCaloSim.cxx.

84{
85
86 #ifdef FCS_DEBUG
87 G4cout<<"[FastCaloSim::ModelTrigger] Got particle with " <<"\n"
88 <<" pdg=" <<fastTrack.GetPrimaryTrack() -> GetDefinition()->GetPDGEncoding() <<"\n"
89 <<" Ekin="<<fastTrack.GetPrimaryTrack() -> GetKineticEnergy() <<"\n"
90 <<" p=" <<fastTrack.GetPrimaryTrack() -> GetMomentum().mag() <<"\n"
91 <<" x=" <<fastTrack.GetPrimaryTrack() -> GetPosition().x() <<"\n"
92 <<" y=" <<fastTrack.GetPrimaryTrack() -> GetPosition().y() <<"\n"
93 <<" z=" <<fastTrack.GetPrimaryTrack() -> GetPosition().z() <<"\n"
94 <<" r=" <<fastTrack.GetPrimaryTrack() -> GetPosition().perp() <<"\n"
95 <<" eta=" <<fastTrack.GetPrimaryTrack() -> GetMomentum().eta() <<"\n"
96 <<" phi=" <<fastTrack.GetPrimaryTrack() -> GetMomentum().phi() <<"\n"
97 <<G4endl;
98 #endif
99
100
101 // Simulate particles below 50 keV with Geant4 to have same config as ISF implementation
102 if (fastTrack.GetPrimaryTrack() -> GetKineticEnergy() < 0.05) {
103 #ifdef FCS_DEBUG
104 G4cout<<"[FastCaloSim::ModelTrigger] Particle below 50 keV threshold. Passing to G4. "<<G4endl;
105 #endif
106 return false;
107 }
108
109 // Check if triggered particle is really on the ID-Calo (parametrization) boundary
110 if (!passedIDCaloBoundary(fastTrack)) {
111 #ifdef FCS_DEBUG
112 G4cout<<"[FastCaloSim::ModelTrigger] Particle failed passedIDCaloBoundary z="<<fastTrack.GetPrimaryTrack() -> GetPosition().z()<<" r="<<fastTrack.GetPrimaryTrack() -> GetPosition().perp()<<G4endl;
113 #endif
114 return false;
115 }
116
117 // Set minimum kinetic energy of pions and other hadrons required to be passed to FastCaloSim
118 float minEkinPions = 200;
119 float minEkinOtherHadrons = 400;
120
121 // Get particle definition
122 const G4ParticleDefinition * G4Particle = fastTrack.GetPrimaryTrack() -> GetDefinition();
123 // Get particle kinetic energy
124 const float Ekin = fastTrack.GetPrimaryTrack() -> GetKineticEnergy();
125
126 // Check particle type
127 bool isPhoton = G4Particle == G4Gamma::Definition();
128 bool isElectron = G4Particle == G4Electron::Definition();
129 bool isPositron = G4Particle == G4Positron::Definition();
130 bool isPionPlus = G4Particle == G4PionPlus::Definition();
131 bool isPionMinus = G4Particle == G4PionMinus::Definition();
132
133 // Pass all photons, electrons and positrons to FastCaloSim
134 if (isPhoton || isElectron || isPositron){
135 #ifdef FCS_DEBUG
136 G4cout<<"[FastCaloSim::ModelTrigger] Model triggered"<<G4endl;
137 #endif
138 return true;
139 }
140
141 // Require minimum kinetic energy of pions needed to be passed to FastCaloSim
142 if (isPionPlus || isPionMinus){
143 bool passMinEkinPions = Ekin > minEkinPions;
144
145 #ifdef FCS_DEBUG
146 if(passMinEkinPions) G4cout<<"[FastCaloSim::ModelTrigger] Model triggered"<<G4endl;
147 else G4cout<<"[FastCaloSim::ModelTrigger] Pion with Ekin="<<Ekin<<" below the minimum "<<minEkinPions<<" MeV threshold. Model not triggered."<<G4endl;
148 #endif
149
150 return passMinEkinPions;
151
152 }
153
154 // Require minimum kinetic energy of other hadrons needed to be passed to FastCaloSim
155 bool passMinEkinOtherHadrons = Ekin > minEkinOtherHadrons;
156
157 #ifdef FCS_DEBUG
158 if(passMinEkinOtherHadrons) G4cout<<"[FastCaloSim::ModelTrigger] Model triggered"<<G4endl;
159 else G4cout<<"[FastCaloSim::ModelTrigger] Other hadron with Ekin="<<Ekin<<" below the minimum "<<minEkinOtherHadrons<<" MeV threshold. Model not triggered."<<G4endl;
160 #endif
161
162 return passMinEkinOtherHadrons;
163}
G4bool passedIDCaloBoundary(const G4FastTrack &fastTrack)
Check if the particle is located at the proper ID-Calo parametrization boundary and is travelling out...

◆ passedIDCaloBoundary()

G4bool FastCaloSim::passedIDCaloBoundary ( const G4FastTrack & fastTrack)

Check if the particle is located at the proper ID-Calo parametrization boundary and is travelling outwards from the ID to the CALO.

Definition at line 296 of file FastCaloSim.cxx.

296 {
297
298
299 /* Method checks if particle has crossed the ID-Calo boundary, defined using three cylinders with pairs of r/z values
300 We also perform a directional check to make sure that the particle does not originate from any backscatter from the MS or CALO */
301
302 // NOTE:
303 // For the inner beam pipe section, innerBeamPipeZ = 4185 corresponds to original AFII boundary, while
304 // innerBeamPipeZ = 4587 corresponds to the CALO::CALO boundary which should be used instead as there
305 // is no triggering volume at this point and we will trigger slightly later than the AFII boundary, so
306 // passedIDCaloBoundary would fail and we would simulate this part with Geant4
307
308 // Barrel
309 const float barrelR = 1148;
310 const float barrelZ = 3549.5;
311 // Inner beam pipe section
312 const float innerBeamPipeR = 120;
313 const float innerBeamPipeZ = 4587;
314 // Outer beam pipe section
315 const float outerBeamPipeR = 41;
316 const float outerBeamPipeZ = 6783;
317
318 // Get particle position
319 const G4ThreeVector particlePosition = fastTrack.GetPrimaryTrack() -> GetPosition();
320 // Get r and z values of position
321 const float r = particlePosition.perp();
322 const float z = particlePosition.z();
323 // Get direction of particle
324 const G4ThreeVector particleDirection = fastTrack.GetPrimaryTrack() -> GetMomentum();
325 // Construct the helper line to decide if particle comes from ID or is backscatter from CALO
326 const G4ThreeVector helperLine = particlePosition + particleDirection;
327
328 // Set 5cm trigger tolerance, i.e. the 'width' of the trigger boundary
329 // be careful not to set this too low, else Geant4 might overjump your trigger boundary and simulate everything with G4
330 const float triggerTolerance = 50;
331
332 // Barrel boundary (horizontal surfaces)
333 if (std::abs(z) <= barrelZ + triggerTolerance) {
334 if (r >= barrelR && r < barrelR + triggerTolerance) return helperLine.perp() >= barrelR;
335 }
336 // Beam pipe boundary (horizontal surfaces)
337 if (std::abs(z) >= barrelZ && std::abs(z) <= innerBeamPipeZ + triggerTolerance){
338 if (r >= innerBeamPipeR && r < innerBeamPipeR + triggerTolerance) return helperLine.perp() >= innerBeamPipeR;
339 }
340 if (std::abs(z) >= innerBeamPipeZ && std::abs(z) <= outerBeamPipeZ){
341 if (r >= outerBeamPipeR && r < outerBeamPipeR + triggerTolerance) return helperLine.perp() >= outerBeamPipeR;
342 }
343 // Barrel boundary (vertical surfaces)
344 if (r >= outerBeamPipeR && r<= innerBeamPipeR){
345 if (std::abs(z) >= innerBeamPipeZ && std::abs(z) < innerBeamPipeZ + triggerTolerance) return std::abs(helperLine.z()) >= innerBeamPipeZ;
346 }
347 // Beam pipe boundary (vertical surfaces)
348 if (r >= innerBeamPipeR && r <= barrelR){
349 if (std::abs(z) >= barrelZ && std::abs(z) < barrelZ + triggerTolerance) return std::abs(helperLine.z()) >= barrelZ;
350 }
351
352 return false;
353
354}
int r
Definition globals.cxx:22

Member Data Documentation

◆ m_CaloCellContainerSDName

std::string FastCaloSim::m_CaloCellContainerSDName
private

Definition at line 75 of file FastCaloSim.h.

◆ m_doG4Transport

bool FastCaloSim::m_doG4Transport
private

Definition at line 77 of file FastCaloSim.h.

◆ m_doPunchThrough

bool FastCaloSim::m_doPunchThrough
private

Definition at line 80 of file FastCaloSim.h.

◆ m_FastCaloSimCaloExtrapolation

PublicToolHandle<IFastCaloSimCaloExtrapolation> FastCaloSim::m_FastCaloSimCaloExtrapolation
private

Definition at line 66 of file FastCaloSim.h.

◆ m_FastCaloSimCaloTransportation

PublicToolHandle<IFastCaloSimCaloTransportation> FastCaloSim::m_FastCaloSimCaloTransportation
private

Definition at line 64 of file FastCaloSim.h.

◆ m_FastCaloSimSvc

ServiceHandle<ISF::IFastCaloSimParamSvc> FastCaloSim::m_FastCaloSimSvc
private

Definition at line 73 of file FastCaloSim.h.

◆ m_G4CaloTransportTool

PublicToolHandle<IG4CaloTransportTool> FastCaloSim::m_G4CaloTransportTool
private

Definition at line 68 of file FastCaloSim.h.

◆ m_PunchThroughSimWrapper

PublicToolHandle<IPunchThroughSimWrapper> FastCaloSim::m_PunchThroughSimWrapper
private

Definition at line 70 of file FastCaloSim.h.


The documentation for this class was generated from the following files: