CMS 3D CMS Logo

CMSHepEmTrackingManager.cc
Go to the documentation of this file.
2 
3 #include "G4EventManager.hh"
4 #include "G4ProcessManager.hh"
5 #include "G4StackManager.hh"
6 #include "G4TrackingManager.hh"
7 
8 CMSHepEmTrackingManager::CMSHepEmTrackingManager(G4double highEnergyLimit) : fHighEnergyLimit(highEnergyLimit) {}
9 
11 
12 void CMSHepEmTrackingManager::BuildPhysicsTable(const G4ParticleDefinition& part) {
13  G4HepEmTrackingManager::BuildPhysicsTable(part);
14 
15  G4ProcessManager* pManager = part.GetProcessManager();
16  G4ProcessManager* pManagerShadow = part.GetMasterProcessManager();
17 
18  G4ProcessVector* pVector = pManager->GetProcessList();
19  for (std::size_t j = 0; j < pVector->size(); ++j) {
20  if (pManagerShadow == pManager) {
21  (*pVector)[j]->BuildPhysicsTable(part);
22  } else {
23  (*pVector)[j]->BuildWorkerPhysicsTable(part);
24  }
25  }
26 }
27 
28 void CMSHepEmTrackingManager::PreparePhysicsTable(const G4ParticleDefinition& part) {
29  G4HepEmTrackingManager::PreparePhysicsTable(part);
30 
31  G4ProcessManager* pManager = part.GetProcessManager();
32  G4ProcessManager* pManagerShadow = part.GetMasterProcessManager();
33 
34  G4ProcessVector* pVector = pManager->GetProcessList();
35  for (std::size_t j = 0; j < pVector->size(); ++j) {
36  if (pManagerShadow == pManager) {
37  (*pVector)[j]->PreparePhysicsTable(part);
38  } else {
39  (*pVector)[j]->PrepareWorkerPhysicsTable(part);
40  }
41  }
42 }
43 
45  if (aTrack->GetKineticEnergy() < fHighEnergyLimit) {
46  // Fully track with G4HepEm.
47  G4HepEmTrackingManager::HandOverOneTrack(aTrack);
48  } else {
49  // Track with the Geant4 kernel and all registered processes.
50  G4EventManager* eventManager = G4EventManager::GetEventManager();
51  G4TrackingManager* trackManager = eventManager->GetTrackingManager();
52 
53  trackManager->ProcessOneTrack(aTrack);
54  if (aTrack->GetTrackStatus() != fStopAndKill) {
55  G4Exception("CMSHepEmTrackingManager::HandOverOneTrack", "NotStopped", FatalException, "track was not stopped");
56  }
57 
58  G4TrackVector* secondaries = trackManager->GimmeSecondaries();
59  eventManager->StackTracks(secondaries);
60  delete aTrack;
61  }
62 }
void PreparePhysicsTable(const G4ParticleDefinition &) override
~CMSHepEmTrackingManager() override
void HandOverOneTrack(G4Track *aTrack) override
part
Definition: HCALResponse.h:20
void BuildPhysicsTable(const G4ParticleDefinition &) override
CMSHepEmTrackingManager(G4double highEnergyLimit)