Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FastSimulationManagerProcess.cc
Go to the documentation of this file.
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//
28//
29//---------------------------------------------------------------
30//
31// G4FastSimulationProcess.cc
32//
33// Description:
34// The process that triggers the parameterised simulations,
35// if any.
36//
37// History:
38// August 97: First implementation. Verderi && MoraDeFreitas.
39// October 06: move to parallel geometry scheme, M. Verderi
40//---------------------------------------------------------------
41
43
46#include "G4ParticleChange.hh"
47#include "G4PathFinder.hh"
49#include "G4ios.hh"
50
52 G4ProcessType theType)
53 : G4VProcess(processName, theType),
54 fWorldVolume(nullptr),
55 fIsTrackingTime(false),
56 fIsFirstStep(false),
57 fGhostNavigator(nullptr),
58 fGhostNavigatorIndex(-1),
59 fIsGhostGeometry(false),
60 fGhostSafety(-1.0),
61 fFieldTrack('0'),
62 fFastSimulationManager(nullptr),
63 fFastSimulationTrigger(false)
64{
65 // -- set Process Sub Type
67
68 fPathFinder = G4PathFinder::GetInstance();
70
71 SetWorldVolume(fTransportationManager->GetNavigatorForTracking()->GetWorldVolume()->GetName());
72 if (verboseLevel > 0)
73 G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
74 << "' is created, and will message geometry with world volume `"
75 << fWorldVolume->GetName() << "'." << G4endl;
77}
78
80 const G4String& worldVolumeName,
81 G4ProcessType theType)
82 : G4VProcess(processName, theType),
83 fWorldVolume(nullptr),
84 fIsTrackingTime(false),
85 fIsFirstStep(false),
86 fGhostNavigator(nullptr),
87 fGhostNavigatorIndex(-1),
88 fIsGhostGeometry(false),
89 fGhostSafety(-1.0),
90 fFieldTrack('0'),
91 fFastSimulationManager(nullptr),
92 fFastSimulationTrigger(false)
93{
94 // -- set Process Sub Type
96
97 fPathFinder = G4PathFinder::GetInstance();
99
100 SetWorldVolume(worldVolumeName);
101 if (verboseLevel > 0)
102 G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
103 << "' is created, and will message geometry with world volume `"
104 << fWorldVolume->GetName() << "'." << G4endl;
106}
107
109 G4VPhysicalVolume* worldVolume,
110 G4ProcessType theType)
111 : G4VProcess(processName, theType),
112 fWorldVolume(nullptr),
113 fIsTrackingTime(false),
114 fIsFirstStep(false),
115 fGhostNavigator(nullptr),
116 fGhostNavigatorIndex(-1),
117 fIsGhostGeometry(false),
118 fGhostSafety(-1.0),
119 fFieldTrack('0'),
120 fFastSimulationManager(nullptr),
121 fFastSimulationTrigger(false)
122{
123 // -- set Process Sub Type
124 SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
125
126 fPathFinder = G4PathFinder::GetInstance();
127 fTransportationManager = G4TransportationManager::GetTransportationManager();
128
129 SetWorldVolume(worldVolume);
130 if (verboseLevel > 0)
131 G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
132 << "' is created, and will message geometry with world volume `"
133 << fWorldVolume->GetName() << "'." << G4endl;
135}
136
141
142// -----------------------
143// User access methods:
144// -----------------------
146{
147 if (fIsTrackingTime) {
149 ed << "G4FastSimulationManagerProcess `" << GetProcessName()
150 << "': changing of world volume at tracking time is not allowed." << G4endl;
151 G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)", "FastSim002",
152 JustWarning, ed, "Call ignored.");
153 }
154 else {
155 G4VPhysicalVolume* newWorld = fTransportationManager->IsWorldExisting(newWorldName);
156 if (newWorld == nullptr) {
157 G4ExceptionDescription tellWhatIsWrong;
158 tellWhatIsWrong << "Volume newWorldName = `" << newWorldName
159 << "' is not a parallel world nor the mass world volume." << G4endl;
160 G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)", "FastSim003",
161 FatalException, tellWhatIsWrong);
162 }
163 if (verboseLevel > 0) {
164 if (fWorldVolume != nullptr)
165 G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
166 << "': changing world volume from '" << fWorldVolume->GetName() << "' to `"
167 << newWorld << "'." << G4endl;
168 else
169 G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
170 << "': setting world volume from to `" << newWorld->GetName() << "'." << G4endl;
171 }
172 fWorldVolume = newWorld;
173 }
174}
175
177{
178 if (newWorld != nullptr)
179 SetWorldVolume(newWorld->GetName());
180 else {
181 G4ExceptionDescription tellWhatIsWrong;
182 tellWhatIsWrong << "Null pointer passed for world volume." << G4endl;
183 G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4VPhysicalVolume* newWorld)",
184 "FastSim004", FatalException, tellWhatIsWrong);
185 }
186}
187
188// --------------------
189// Start/End tracking:
190// --------------------
192{
193 fIsTrackingTime = true;
194 fIsFirstStep = true;
195
196 // -- fetch the navigator (and its index) and activate it:
197 G4TransportationManager* transportationManager =
199 fGhostNavigator = transportationManager->GetNavigator(fWorldVolume);
200 fIsGhostGeometry = (fGhostNavigator != transportationManager->GetNavigatorForTracking());
201 if (fIsGhostGeometry)
202 fGhostNavigatorIndex = transportationManager->ActivateNavigator(fGhostNavigator);
203 else
204 fGhostNavigatorIndex = -1;
205
206 fPathFinder->PrepareNewTrack(track->GetPosition(), track->GetMomentumDirection());
207}
208
210{
211 fIsTrackingTime = false;
212 if (fIsGhostGeometry) fTransportationManager->DeActivateNavigator(fGhostNavigator);
213}
214
215// ------------------------------------------
216// PostStepGetPhysicalInteractionLength():
217// ------------------------------------------
221{
222 // -- Get current volume, and check for presence of fast simulation manager.
223 // -- For the case of the navigator for tracking (fGhostNavigatorIndex == 0)
224 // -- we use the track volume. This allows the code to be valid for both
225 // -- cases where the PathFinder is used (G4CoupledTranportation) or not
226 // -- (G4Transportation).
227 const G4VPhysicalVolume* currentVolume(nullptr);
228 if (fIsGhostGeometry)
229 currentVolume = fPathFinder->GetLocatedVolume(fGhostNavigatorIndex);
230 else
231 currentVolume = track.GetVolume();
232
233 if (currentVolume != nullptr) {
234 fFastSimulationManager = currentVolume->GetLogicalVolume()->GetFastSimulationManager();
235 if (fFastSimulationManager != nullptr) {
236 // Ask for trigger:
237 fFastSimulationTrigger =
238 fFastSimulationManager->PostStepGetFastSimulationManagerTrigger(track, fGhostNavigator);
239 if (fFastSimulationTrigger) {
240 // Take control over stepping:
242 return 0.0;
243 }
244 }
245 }
246
247 // -- no fast simulation occuring there:
249 return DBL_MAX;
250}
251
252//------------------------------------
253// PostStepDoIt()
254//------------------------------------
256{
257 G4VParticleChange* finalState = fFastSimulationManager->InvokePostStepDoIt();
258
259 // If the particle is still alive, suspend it to force physics re-initialisation:
260 if (finalState->GetTrackStatus() != fStopAndKill) finalState->ProposeTrackStatus(fSuspend);
261
262 return finalState;
263}
264
266 const G4Track& track, G4double previousStepSize, G4double currentMinimumStep,
267 G4double& proposedSafety, G4GPILSelection* selection)
268{
269 *selection = NotCandidateForSelection;
270 G4double returnedStep = DBL_MAX;
271
272 // ---------------------------------------------------
273 // -- Below code valid for ghost geometry, otherwise
274 // -- useless for fast simulation attached to mass
275 // -- geometry. Warn user in case along used for
276 // -- mass geometry ?
277 // --------------------------------------------------
278 if (fIsGhostGeometry) {
279 static G4ThreadLocal G4FieldTrack* endTrack_G4MT_TLS_ = nullptr;
280 if (endTrack_G4MT_TLS_ == nullptr) endTrack_G4MT_TLS_ = new G4FieldTrack('0');
281 G4FieldTrack& endTrack = *endTrack_G4MT_TLS_;
282
283 static G4ThreadLocal ELimited* eLimited_G4MT_TLS_ = nullptr;
284 if (eLimited_G4MT_TLS_ == nullptr) eLimited_G4MT_TLS_ = new ELimited;
285 ELimited& eLimited = *eLimited_G4MT_TLS_;
286
287 if (previousStepSize > 0.) fGhostSafety -= previousStepSize;
288 if (fGhostSafety < 0.) fGhostSafety = 0.0;
289
290 // ------------------------------------------
291 // Determination of the proposed step length:
292 // ------------------------------------------
293 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.) {
294 // -- No chance to limit the step, as proposed move inside safety
295 returnedStep = currentMinimumStep;
296 proposedSafety = fGhostSafety - currentMinimumStep;
297 }
298 else {
299 // -- Proposed move exceeds safety, need to state
300 G4FieldTrackUpdator::Update(&fFieldTrack, &track);
301 returnedStep = fPathFinder->ComputeStep(fFieldTrack, currentMinimumStep, fGhostNavigatorIndex,
302 track.GetCurrentStepNumber(), fGhostSafety, eLimited,
303 endTrack, track.GetVolume());
304
305 if (eLimited == kDoNot)
306 fGhostSafety =
307 fGhostNavigator->ComputeSafety(endTrack.GetPosition()); // -- step no limited by ghost
308 proposedSafety = fGhostSafety;
309 if (eLimited == kUnique || eLimited == kSharedOther)
310 *selection = CandidateForSelection;
311 else if (eLimited == kSharedTransport)
312 returnedStep *=
313 (1.0 + 1.0e-9); // -- Expand to disable its selection in Step Manager comparison
314 }
315 }
316
317 // ----------------------------------------------
318 // Returns the fGhostSafety as the proposedSafety
319 // The SteppingManager will take care of keeping
320 // the smallest one.
321 // ----------------------------------------------
322 return returnedStep;
323}
324
326 const G4Step&)
327{
328 fDummyParticleChange.Initialize(track);
329 return &fDummyParticleChange;
330}
331
332//--------------------------------------------
333// At Rest parameterisation:
334//--------------------------------------------
335// AtRestGetPhysiscalInteractionLength:
336//--------------------------------------------
340{
341 const G4VPhysicalVolume* currentVolume(nullptr);
342 if (fIsGhostGeometry)
343 currentVolume = fPathFinder->GetLocatedVolume(fGhostNavigatorIndex);
344 else
345 currentVolume = track.GetVolume();
346 fFastSimulationManager = currentVolume->GetLogicalVolume()->GetFastSimulationManager();
347 if (fFastSimulationManager != nullptr) {
348 // Ask for trigger:
349 fFastSimulationTrigger =
350 fFastSimulationManager->AtRestGetFastSimulationManagerTrigger(track, fGhostNavigator);
351 if (fFastSimulationTrigger) {
352 // Dirty trick to take control over stepping. Does anyone will ever use that ?
354 return -1.0;
355 }
356 }
357
358 // -- no fast simulation occuring there:
360 return DBL_MAX;
361}
362
363//-----------------------------------------------
364// AtRestDoIt:
365//-----------------------------------------------
367{
368 return fFastSimulationManager->InvokeAtRestDoIt();
369}
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4ForceCondition
@ NotForced
@ ExclusivelyForced
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
@ kDoNot
@ kUnique
@ kSharedOther
@ kSharedTransport
G4ProcessType
@ fSuspend
@ fStopAndKill
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *) override
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4FastSimulationManagerProcess(const G4String &processName="G4FastSimulationManagerProcess", G4ProcessType theType=fParameterisation)
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &) override
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step) override
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection) override
G4VParticleChange * InvokePostStepDoIt()
G4VParticleChange * InvokeAtRestDoIt()
G4bool PostStepGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=nullptr)
G4bool AtRestGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=nullptr)
static void Update(G4FieldTrack *, const G4Track *)
G4ThreeVector GetPosition() const
void RemoveFSMP(G4FastSimulationManagerProcess *)
void AddFSMP(G4FastSimulationManagerProcess *)
static G4GlobalFastSimulationManager * GetGlobalFastSimulationManager()
G4FastSimulationManager * GetFastSimulationManager() const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4VPhysicalVolume * GetWorldVolume() const
G4VPhysicalVolume * GetLocatedVolume(G4int navId) const
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=nullptr)
static G4PathFinder * GetInstance()
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4int GetCurrentStepNumber() const
const G4ThreeVector & GetMomentumDirection() const
static G4TransportationManager * GetTransportationManager()
G4VPhysicalVolume * IsWorldExisting(const G4String &worldName)
G4Navigator * GetNavigatorForTracking() const
G4int ActivateNavigator(G4Navigator *aNavigator)
void DeActivateNavigator(G4Navigator *aNavigator)
G4Navigator * GetNavigator(const G4String &worldName)
void ProposeTrackStatus(G4TrackStatus status)
virtual void Initialize(const G4Track &)
G4TrackStatus GetTrackStatus() const
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const
#define DBL_MAX
Definition templates.hh:62
#define G4ThreadLocal
Definition tls.hh:77