Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/runAndEvent/RE07/src/SpecializedTrackingManager.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 //
 27 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 28 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 29 
 30 #include "SpecializedTrackingManager.hh"
 31 
 32 #include "G4EventManager.hh"
 33 #include "G4ProcessManager.hh"
 34 #include "G4RegionStore.hh"
 35 #include "G4StackManager.hh"
 36 #include "G4SystemOfUnits.hh"
 37 #include "G4TrackingManager.hh"
 38 
 39 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 40 
 41 SpecializedTrackingManager::SpecializedTrackingManager() {}
 42 
 43 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 44 
 45 SpecializedTrackingManager::~SpecializedTrackingManager() {}
 46 
 47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 48 
 49 void SpecializedTrackingManager::BuildPhysicsTable(const G4ParticleDefinition& part)
 50 {
 51   if (fBackRegion == nullptr) {
 52     fBackRegion = G4RegionStore::GetInstance()->GetRegion("Back", false);
 53   }
 54 
 55   G4ProcessManager* pManager = part.GetProcessManager();
 56   G4ProcessManager* pManagerShadow = part.GetMasterProcessManager();
 57 
 58   G4ProcessVector* pVector = pManager->GetProcessList();
 59   for (std::size_t j = 0; j < pVector->size(); ++j) {
 60     if (pManagerShadow == pManager) {
 61       (*pVector)[j]->BuildPhysicsTable(part);
 62     }
 63     else {
 64       (*pVector)[j]->BuildWorkerPhysicsTable(part);
 65     }
 66   }
 67 }
 68 
 69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 70 
 71 void SpecializedTrackingManager::PreparePhysicsTable(const G4ParticleDefinition& part)
 72 {
 73   G4ProcessManager* pManager = part.GetProcessManager();
 74   G4ProcessManager* pManagerShadow = part.GetMasterProcessManager();
 75 
 76   G4ProcessVector* pVector = pManager->GetProcessList();
 77   for (std::size_t j = 0; j < pVector->size(); ++j) {
 78     if (pManagerShadow == pManager) {
 79       (*pVector)[j]->PreparePhysicsTable(part);
 80     }
 81     else {
 82       (*pVector)[j]->PrepareWorkerPhysicsTable(part);
 83     }
 84   }
 85 }
 86 
 87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 88 
 89 void SpecializedTrackingManager::HandOverOneTrack(G4Track* aTrack)
 90 {
 91   if (aTrack->GetKineticEnergy() < 100 * MeV) {
 92     // If the particle energy is lower than 100 MeV, track it immediately by
 93     // passing to the generic G4TrackingManager. This avoids storing lower
 94     // energy particles in the buffer and feeding it through the specialized
 95     // tracking.
 96     G4EventManager* eventManager = G4EventManager::GetEventManager();
 97     G4TrackingManager* trackManager = eventManager->GetTrackingManager();
 98 
 99     trackManager->ProcessOneTrack(aTrack);
100     if (aTrack->GetTrackStatus() != fStopAndKill) {
101       G4Exception("SpecializedTrackingManager::HandOverOneTrack", "NotStopped", FatalException,
102                   "track was not stopped");
103     }
104 
105     G4TrackVector* secondaries = trackManager->GimmeSecondaries();
106     eventManager->StackTracks(secondaries);
107     delete aTrack;
108     return;
109   }
110 
111   fBufferedTracks.push_back(aTrack);
112 }
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
115 
116 void SpecializedTrackingManager::FlushEvent()
117 {
118   G4EventManager* eventManager = G4EventManager::GetEventManager();
119   G4TrackingManager* trackManager = eventManager->GetTrackingManager();
120   G4SteppingManager* steppingManager = trackManager->GetSteppingManager();
121   G4TrackVector* secondaries = trackManager->GimmeSecondaries();
122 
123   for (G4Track* aTrack : fBufferedTracks) {
124     // Clear secondary particle vector
125     for (std::size_t itr = 0; itr < secondaries->size(); ++itr) {
126       delete (*secondaries)[itr];
127     }
128     secondaries->clear();
129 
130     steppingManager->SetInitialStep(aTrack);
131 
132     G4UserTrackingAction* userTrackingAction = trackManager->GetUserTrackingAction();
133     if (userTrackingAction != nullptr) {
134       userTrackingAction->PreUserTrackingAction(aTrack);
135     }
136 
137     // Give SteppingManger the maxmimum number of processes
138     steppingManager->GetProcessNumber();
139 
140     // Give track the pointer to the Step
141     aTrack->SetStep(steppingManager->GetStep());
142 
143     // Inform beginning of tracking to physics processes
144     aTrack->GetDefinition()->GetProcessManager()->StartTracking(aTrack);
145 
146     // Track the particle Step-by-Step while it is alive
147     while ((aTrack->GetTrackStatus() == fAlive) || (aTrack->GetTrackStatus() == fStopButAlive)) {
148       G4Region* region = aTrack->GetVolume()->GetLogicalVolume()->GetRegion();
149       if (region == fBackRegion) {
150         StepInBackRegion(aTrack);
151       }
152       else {
153         StepOutside(aTrack);
154       }
155     }
156 
157     aTrack->GetDefinition()->GetProcessManager()->EndTracking();
158 
159     if (userTrackingAction != nullptr) {
160       userTrackingAction->PostUserTrackingAction(aTrack);
161     }
162 
163     eventManager->StackTracks(secondaries);
164     delete aTrack;
165   }
166 
167   fBufferedTracks.clear();
168 }
169 
170 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
171 
172 void SpecializedTrackingManager::StepInBackRegion(G4Track* aTrack)
173 {
174   G4EventManager* eventManager = G4EventManager::GetEventManager();
175   G4TrackingManager* trackManager = eventManager->GetTrackingManager();
176   G4SteppingManager* steppingManager = trackManager->GetSteppingManager();
177 
178   // Track the particle Step-by-Step while it is alive and inside the "Back"
179   // region of the detector. Implement a low-energy cut-off for particles
180   // below 100 MeV. More specialized handling would also be possible, such
181   // as only killing particles in non-sensitive materials / volumes.
182   while ((aTrack->GetTrackStatus() == fAlive) || (aTrack->GetTrackStatus() == fStopButAlive)) {
183     aTrack->IncrementCurrentStepNumber();
184     steppingManager->Stepping();
185 
186     if (aTrack->GetTrackStatus() != fStopAndKill) {
187       // Switch the touchable to update the volume, which is checked in the
188       // condition below and at the call site.
189       aTrack->SetTouchableHandle(aTrack->GetNextTouchableHandle());
190       G4Region* region = aTrack->GetVolume()->GetLogicalVolume()->GetRegion();
191       if (region != fBackRegion) {
192         return;
193       }
194 
195       if (aTrack->GetKineticEnergy() < 100 * MeV) {
196         // Kill the particle.
197         aTrack->SetTrackStatus(fStopAndKill);
198       }
199     }
200   }
201 }
202 
203 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
204 
205 void SpecializedTrackingManager::StepOutside(G4Track* aTrack)
206 {
207   G4EventManager* eventManager = G4EventManager::GetEventManager();
208   G4TrackingManager* trackManager = eventManager->GetTrackingManager();
209   G4SteppingManager* steppingManager = trackManager->GetSteppingManager();
210 
211   // Track the particle Step-by-Step while it is alive and still outside of
212   // the "Back" region.
213   while ((aTrack->GetTrackStatus() == fAlive) || (aTrack->GetTrackStatus() == fStopButAlive)) {
214     aTrack->IncrementCurrentStepNumber();
215     steppingManager->Stepping();
216 
217     if (aTrack->GetTrackStatus() != fStopAndKill) {
218       // Switch the touchable to update the volume, which is checked in the
219       // condition below and at the call site.
220       aTrack->SetTouchableHandle(aTrack->GetNextTouchableHandle());
221       G4Region* region = aTrack->GetVolume()->GetLogicalVolume()->GetRegion();
222       if (region == fBackRegion) {
223         return;
224       }
225     }
226   }
227 }
228