Geant4 11.1.1
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
42#include "G4ios.hh"
46#include "G4PathFinder.hh"
47#include "G4ParticleChange.hh"
49
50#define PARANOIA
51
54 G4ProcessType theType) :
55 G4VProcess(processName,theType),
56 fWorldVolume ( nullptr ),
57 fIsTrackingTime ( false ),
58 fIsFirstStep ( false ),
59 fGhostNavigator ( nullptr ),
60 fGhostNavigatorIndex ( -1 ),
61 fIsGhostGeometry ( false ),
62 fGhostSafety ( -1.0 ),
63 fFieldTrack ( '0' ),
64 fFastSimulationManager( nullptr ),
65 fFastSimulationTrigger( false )
66{
67 // -- set Process Sub Type
69
70
71 fPathFinder = G4PathFinder::GetInstance();
73
74 SetWorldVolume(fTransportationManager->GetNavigatorForTracking()->GetWorldVolume()->GetName());
75 if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
76 << "' is created, and will message geometry with world volume `"
77 << fWorldVolume->GetName() << "'." << G4endl;
79}
80
81
84 const G4String& worldVolumeName,
85 G4ProcessType theType) :
86 G4VProcess(processName,theType),
87 fWorldVolume ( nullptr ),
88 fIsTrackingTime ( false ),
89 fIsFirstStep ( false ),
90 fGhostNavigator ( nullptr ),
91 fGhostNavigatorIndex ( -1 ),
92 fIsGhostGeometry ( false ),
93 fGhostSafety ( -1.0 ),
94 fFieldTrack ( '0' ),
95 fFastSimulationManager( nullptr ),
96 fFastSimulationTrigger( false )
97{
98 // -- set Process Sub Type
100
101
102 fPathFinder = G4PathFinder::GetInstance();
103 fTransportationManager = G4TransportationManager::GetTransportationManager();
104
105 SetWorldVolume(worldVolumeName);
106 if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
107 << "' is created, and will message geometry with world volume `"
108 << fWorldVolume->GetName() << "'." << G4endl;
110}
111
112
115 G4VPhysicalVolume* worldVolume,
116 G4ProcessType theType) :
117 G4VProcess(processName,theType),
118 fWorldVolume ( nullptr ),
119 fIsTrackingTime ( false ),
120 fIsFirstStep ( false ),
121 fGhostNavigator ( nullptr ),
122 fGhostNavigatorIndex ( -1 ),
123 fIsGhostGeometry ( false ),
124 fGhostSafety ( -1.0 ),
125 fFieldTrack ( '0' ),
126 fFastSimulationManager( nullptr ),
127 fFastSimulationTrigger( false )
128{
129 // -- set Process Sub Type
130 SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
131
132
133 fPathFinder = G4PathFinder::GetInstance();
134 fTransportationManager = G4TransportationManager::GetTransportationManager();
135
136 SetWorldVolume(worldVolume);
137 if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
138 << "' is created, and will message geometry with world volume `"
139 << fWorldVolume->GetName() << "'." << G4endl;
141}
142
143
145{
147}
148
149
150// -----------------------
151// User access methods:
152// -----------------------
154{
155 if (fIsTrackingTime)
156 {
158 ed << "G4FastSimulationManagerProcess `" << GetProcessName()
159 << "': changing of world volume at tracking time is not allowed." << G4endl;
160 G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
161 "FastSim002",
162 JustWarning, ed,
163 "Call ignored.");
164 }
165 else
166 {
167 G4VPhysicalVolume* newWorld = fTransportationManager->IsWorldExisting(newWorldName);
168 if (newWorld == 0)
169 {
170 G4ExceptionDescription tellWhatIsWrong;
171 tellWhatIsWrong << "Volume newWorldName = `" << newWorldName
172 << "' is not a parallel world nor the mass world volume."
173 << G4endl;
174 G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
175 "FastSim003",
177 tellWhatIsWrong);
178 }
179 if (verboseLevel>0)
180 {
181 if (fWorldVolume) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
182 << "': changing world volume from '" << fWorldVolume->GetName()
183 << "' to `" << newWorld << "'." << G4endl;
184 else G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
185 << "': setting world volume from to `"<< newWorld->GetName() << "'." << G4endl;
186 }
187 fWorldVolume = newWorld;
188 }
189}
190
191
193{
194 if (newWorld) SetWorldVolume(newWorld->GetName());
195 else
196 {
197 G4ExceptionDescription tellWhatIsWrong;
198 tellWhatIsWrong << "Null pointer passed for world volume." << G4endl;
199 G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4VPhysicalVolume* newWorld)",
200 "FastSim004",
202 tellWhatIsWrong);
203 }
204}
205
206
207// --------------------
208// Start/End tracking:
209// --------------------
211{
212 fIsTrackingTime = true;
213 fIsFirstStep = true;
214
215 // -- fetch the navigator (and its index) and activate it:
217 fGhostNavigator = transportationManager->GetNavigator(fWorldVolume);
218 fIsGhostGeometry = (fGhostNavigator != transportationManager->GetNavigatorForTracking());
219 if (fIsGhostGeometry) fGhostNavigatorIndex = transportationManager->ActivateNavigator(fGhostNavigator);
220 else fGhostNavigatorIndex = -1;
221
222 fPathFinder->PrepareNewTrack(track->GetPosition(), track->GetMomentumDirection());
223}
224
225
226void
229{
230 fIsTrackingTime = false;
231 if ( fIsGhostGeometry ) fTransportationManager->DeActivateNavigator(fGhostNavigator);
232}
233
234
235// ------------------------------------------
236// PostStepGetPhysicalInteractionLength():
237// ------------------------------------------
241 G4double,
243{
244 // -- Get current volume, and check for presence of fast simulation manager.
245 // -- For the case of the navigator for tracking (fGhostNavigatorIndex == 0)
246 // -- we use the track volume. This allows the code to be valid for both
247 // -- cases where the PathFinder is used (G4CoupledTranportation) or not
248 // -- (G4Transportation).
249 const G4VPhysicalVolume* currentVolume(0);
250 if ( fIsGhostGeometry ) currentVolume = fPathFinder->GetLocatedVolume(fGhostNavigatorIndex);
251 else currentVolume = track.GetVolume();
252
253 if ( currentVolume )
254 {
255 fFastSimulationManager = currentVolume->GetLogicalVolume()->GetFastSimulationManager();
256 if( fFastSimulationManager )
257 {
258 // Ask for trigger:
259 fFastSimulationTrigger = fFastSimulationManager->PostStepGetFastSimulationManagerTrigger(track, fGhostNavigator);
260 if( fFastSimulationTrigger )
261 {
262 // Take control over stepping:
264 return 0.0;
265 }
266 }
267 }
268
269 // -- no fast simulation occuring there:
271 return DBL_MAX;
272}
273
274//------------------------------------
275// PostStepDoIt()
276//------------------------------------
279PostStepDoIt(const G4Track&,
280 const G4Step&)
281{
282 G4VParticleChange* finalState = fFastSimulationManager->InvokePostStepDoIt();
283
284 // If the particle is still alive, suspend it to force physics re-initialisation:
285 if (finalState->GetTrackStatus() != fStopAndKill) finalState->ProposeTrackStatus(fSuspend);
286
287 return finalState;
288}
289
290
294 G4double previousStepSize,
295 G4double currentMinimumStep,
296 G4double& proposedSafety,
297 G4GPILSelection* selection)
298{
299
300 *selection = NotCandidateForSelection;
301 G4double returnedStep = DBL_MAX;
302
303 // ---------------------------------------------------
304 // -- Below code valid for ghost geometry, otherwise
305 // -- useless for fast simulation attached to mass
306 // -- geometry. Warn user in case along used for
307 // -- mass geometry ?
308 // --------------------------------------------------
309 if ( fIsGhostGeometry )
310 {
311 static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ;
312 if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ;
313 G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
314
315 static G4ThreadLocal ELimited *eLimited_G4MT_TLS_ = 0 ;
316 if (!eLimited_G4MT_TLS_) eLimited_G4MT_TLS_ = new ELimited ;
317 ELimited &eLimited = *eLimited_G4MT_TLS_;
318
319 if (previousStepSize > 0.) fGhostSafety -= previousStepSize;
320 if (fGhostSafety < 0.) fGhostSafety = 0.0;
321
322 // ------------------------------------------
323 // Determination of the proposed step length:
324 // ------------------------------------------
325 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
326 {
327 // -- No chance to limit the step, as proposed move inside safety
328 returnedStep = currentMinimumStep;
329 proposedSafety = fGhostSafety - currentMinimumStep;
330 }
331 else
332 {
333 // -- Proposed move exceeds safety, need to state
334 G4FieldTrackUpdator::Update(&fFieldTrack, &track);
335 returnedStep = fPathFinder->ComputeStep(fFieldTrack,
336 currentMinimumStep,
337 fGhostNavigatorIndex,
338 track.GetCurrentStepNumber(),
339 fGhostSafety,
340 eLimited,
341 endTrack,
342 track.GetVolume());
343
344 if(eLimited == kDoNot) fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition()); // -- step no limited by ghost
345 proposedSafety = fGhostSafety;
346 if (eLimited == kUnique || eLimited == kSharedOther) *selection = CandidateForSelection;
347 else if (eLimited == kSharedTransport) returnedStep *= (1.0 + 1.0e-9); // -- Expand to disable its selection in Step Manager comparison
348 }
349 }
350
351
352 // ----------------------------------------------
353 // Returns the fGhostSafety as the proposedSafety
354 // The SteppingManager will take care of keeping
355 // the smallest one.
356 // ----------------------------------------------
357 return returnedStep;
358}
359
362AlongStepDoIt(const G4Track& track,
363 const G4Step&)
364{
365 fDummyParticleChange.Initialize(track);
366 return &fDummyParticleChange;
367}
368
369
370//--------------------------------------------
371// At Rest parameterisation:
372//--------------------------------------------
373// AtRestGetPhysiscalInteractionLength:
374//--------------------------------------------
379{
380 const G4VPhysicalVolume* currentVolume(0);
381 if ( fIsGhostGeometry ) currentVolume = fPathFinder->GetLocatedVolume(fGhostNavigatorIndex);
382 else currentVolume = track.GetVolume();
383 fFastSimulationManager = currentVolume->GetLogicalVolume()->GetFastSimulationManager();
384 if( fFastSimulationManager )
385 {
386 // Ask for trigger:
387 fFastSimulationTrigger = fFastSimulationManager->AtRestGetFastSimulationManagerTrigger(track, fGhostNavigator);
388 if( fFastSimulationTrigger )
389 {
390 // Dirty trick to take control over stepping. Does anyone will ever use that ?
392 return -1.0;
393 }
394 }
395
396 // -- no fast simulation occuring there:
398 return DBL_MAX;
399}
400
401//-----------------------------------------------
402// AtRestDoIt:
403//-----------------------------------------------
405{
406 return fFastSimulationManager->InvokeAtRestDoIt();
407}
408
409
411{
412 /* G4cout << " >>>>> Trigger Status : ";
413 switch(fFastSimulationManager->GetTriggerStatus())
414 {
415 case NoModel:
416 G4cout << "NoModel" << G4endl;
417 break;
418 case OnBoundaryButLeaving:
419 G4cout << "OnBoundaryButLeaving" << G4endl;
420 break;
421 case OneModelTrigger:
422 G4cout << "OneModelTrigger" << G4endl;
423 break;
424 case NoModelTrigger:
425 G4cout << "NoModelTrigger" << G4endl;
426 break;
427 case Undefined:
428 G4cout << "Undefined" << G4endl;
429 break;
430 default:
431 G4cout << " Bizarre..." << G4endl;
432 break;
433 }*/
434}
435
436
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ NotForced
@ ExclusivelyForced
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
ELimited
@ kDoNot
@ kUnique
@ kSharedOther
@ kSharedTransport
G4ProcessType
@ fSuspend
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
G4FastSimulationManagerProcess(const G4String &processName="G4FastSimulationManagerProcess", G4ProcessType theType=fParameterisation)
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &step)
G4VParticleChange * InvokePostStepDoIt()
G4VParticleChange * InvokeAtRestDoIt()
G4bool PostStepGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=0)
G4bool AtRestGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=0)
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()
Definition: G4PathFinder.cc:52
Definition: G4Step.hh:62
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
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
#define DBL_MAX
Definition: templates.hh:62
#define G4ThreadLocal
Definition: tls.hh:77