Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4FastSimulationManagerProcess Class Reference

#include <G4FastSimulationManagerProcess.hh>

Inheritance diagram for G4FastSimulationManagerProcess:
Collaboration diagram for G4FastSimulationManagerProcess:

Public Member Functions

 G4FastSimulationManagerProcess (const G4String &processName="G4FastSimulationManagerProcess", G4ProcessType theType=fParameterisation)
 
 G4FastSimulationManagerProcess (const G4String &processName, const G4String &worldVolumeName, G4ProcessType theType=fParameterisation)
 
 G4FastSimulationManagerProcess (const G4String &processName, G4VPhysicalVolume *worldVolume, G4ProcessType theType=fParameterisation)
 
virtual ~G4FastSimulationManagerProcess ()
 
G4VPhysicalVolumeGetWorldVolume () const
 
void SetWorldVolume (G4String)
 
void SetWorldVolume (G4VPhysicalVolume *)
 
void StartTracking (G4Track *)
 
void EndTracking ()
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &step)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
void Verbose () const
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 78 of file G4FastSimulationManagerProcess.hh.

Constructor & Destructor Documentation

G4FastSimulationManagerProcess::G4FastSimulationManagerProcess ( const G4String processName = "G4FastSimulationManagerProcess",
G4ProcessType  theType = fParameterisation 
)

Definition at line 54 of file G4FastSimulationManagerProcess.cc.

55  :
56  G4VProcess(processName,theType),
57  fWorldVolume ( nullptr ),
58  fIsTrackingTime ( false ),
59  fIsFirstStep ( false ),
60  fGhostNavigator ( nullptr ),
61  fGhostNavigatorIndex ( -1 ),
62  fIsGhostGeometry ( false ),
63  fGhostSafety ( -1.0 ),
64  fFieldTrack ( '0' ),
65  fFastSimulationManager( nullptr ),
66  fFastSimulationTrigger( false )
67 {
68  // -- set Process Sub Type
69  SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
70 
71 
72  fPathFinder = G4PathFinder::GetInstance();
73  fTransportationManager = G4TransportationManager::GetTransportationManager();
74 
75  SetWorldVolume(fTransportationManager->GetNavigatorForTracking()->GetWorldVolume()->GetName());
76  if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
77  << "' is created, and will message geometry with world volume `"
78  << fWorldVolume->GetName() << "'." << G4endl;
80 }
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:57
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:52
G4Navigator * GetNavigatorForTracking() const
static G4GlobalFastSimulationManager * GetGlobalFastSimulationManager()
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
static G4TransportationManager * GetTransportationManager()
#define G4endl
Definition: G4ios.hh:61
G4VPhysicalVolume * GetWorldVolume() const
void AddFSMP(G4FastSimulationManagerProcess *)

Here is the call graph for this function:

G4FastSimulationManagerProcess::G4FastSimulationManagerProcess ( const G4String processName,
const G4String worldVolumeName,
G4ProcessType  theType = fParameterisation 
)

Definition at line 84 of file G4FastSimulationManagerProcess.cc.

86  :
87  G4VProcess(processName,theType),
88  fWorldVolume ( nullptr ),
89  fIsTrackingTime ( false ),
90  fIsFirstStep ( false ),
91  fGhostNavigator ( nullptr ),
92  fGhostNavigatorIndex ( -1 ),
93  fIsGhostGeometry ( false ),
94  fGhostSafety ( -1.0 ),
95  fFieldTrack ( '0' ),
96  fFastSimulationManager( nullptr ),
97  fFastSimulationTrigger( false )
98 {
99  // -- set Process Sub Type
100  SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
101 
102 
103  fPathFinder = G4PathFinder::GetInstance();
104  fTransportationManager = G4TransportationManager::GetTransportationManager();
105 
106  SetWorldVolume(worldVolumeName);
107  if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
108  << "' is created, and will message geometry with world volume `"
109  << fWorldVolume->GetName() << "'." << G4endl;
111 }
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:57
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:52
static G4GlobalFastSimulationManager * GetGlobalFastSimulationManager()
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
static G4TransportationManager * GetTransportationManager()
#define G4endl
Definition: G4ios.hh:61
void AddFSMP(G4FastSimulationManagerProcess *)

Here is the call graph for this function:

G4FastSimulationManagerProcess::G4FastSimulationManagerProcess ( const G4String processName,
G4VPhysicalVolume worldVolume,
G4ProcessType  theType = fParameterisation 
)

Definition at line 115 of file G4FastSimulationManagerProcess.cc.

117  :
118  G4VProcess(processName,theType),
119  fWorldVolume ( nullptr ),
120  fIsTrackingTime ( false ),
121  fIsFirstStep ( false ),
122  fGhostNavigator ( nullptr ),
123  fGhostNavigatorIndex ( -1 ),
124  fIsGhostGeometry ( false ),
125  fGhostSafety ( -1.0 ),
126  fFieldTrack ( '0' ),
127  fFastSimulationManager( nullptr ),
128  fFastSimulationTrigger( false )
129 {
130  // -- set Process Sub Type
131  SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
132 
133 
134  fPathFinder = G4PathFinder::GetInstance();
135  fTransportationManager = G4TransportationManager::GetTransportationManager();
136 
137  SetWorldVolume(worldVolume);
138  if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
139  << "' is created, and will message geometry with world volume `"
140  << fWorldVolume->GetName() << "'." << G4endl;
142 }
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:57
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:52
static G4GlobalFastSimulationManager * GetGlobalFastSimulationManager()
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
static G4TransportationManager * GetTransportationManager()
#define G4endl
Definition: G4ios.hh:61
void AddFSMP(G4FastSimulationManagerProcess *)

Here is the call graph for this function:

G4FastSimulationManagerProcess::~G4FastSimulationManagerProcess ( )
virtual

Definition at line 145 of file G4FastSimulationManagerProcess.cc.

146 {
148 }
void RemoveFSMP(G4FastSimulationManagerProcess *)
static G4GlobalFastSimulationManager * GetGlobalFastSimulationManager()

Here is the call graph for this function:

Member Function Documentation

G4VParticleChange * G4FastSimulationManagerProcess::AlongStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4VProcess.

Definition at line 363 of file G4FastSimulationManagerProcess.cc.

365 {
366  fDummyParticleChange.Initialize(track);
367  return &fDummyParticleChange;
368 }
virtual void Initialize(const G4Track &)

Here is the call graph for this function:

G4double G4FastSimulationManagerProcess::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
)
virtual

Implements G4VProcess.

Definition at line 294 of file G4FastSimulationManagerProcess.cc.

299 {
300 
301  *selection = NotCandidateForSelection;
302  G4double returnedStep = DBL_MAX;
303 
304  // ---------------------------------------------------
305  // -- Below code valid for ghost geometry, otherwise
306  // -- useless for fast simulation attached to mass
307  // -- geometry. Warn user in case along used for
308  // -- mass geometry ?
309  // --------------------------------------------------
310  if ( fIsGhostGeometry )
311  {
312  static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ;
313  if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ;
314  G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
315 
316  static G4ThreadLocal ELimited *eLimited_G4MT_TLS_ = 0 ;
317  if (!eLimited_G4MT_TLS_) eLimited_G4MT_TLS_ = new ELimited ;
318  ELimited &eLimited = *eLimited_G4MT_TLS_;
319 
320  if (previousStepSize > 0.) fGhostSafety -= previousStepSize;
321  if (fGhostSafety < 0.) fGhostSafety = 0.0;
322 
323  // ------------------------------------------
324  // Determination of the proposed step length:
325  // ------------------------------------------
326  if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
327  {
328  // -- No chance to limit the step, as proposed move inside safety
329  returnedStep = currentMinimumStep;
330  proposedSafety = fGhostSafety - currentMinimumStep;
331  }
332  else
333  {
334  // -- Proposed move exceeds safety, need to state
335  G4FieldTrackUpdator::Update(&fFieldTrack, &track);
336  returnedStep = fPathFinder->ComputeStep(fFieldTrack,
337  currentMinimumStep,
338  fGhostNavigatorIndex,
339  track.GetCurrentStepNumber(),
340  fGhostSafety,
341  eLimited,
342  endTrack,
343  track.GetVolume());
344 
345  if(eLimited == kDoNot) fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition()); // -- step no limited by ghost
346  proposedSafety = fGhostSafety;
347  if (eLimited == kUnique || eLimited == kSharedOther) *selection = CandidateForSelection;
348  else if (eLimited == kSharedTransport) returnedStep *= (1.0 + 1.0e-9); // -- Expand to disable its selection in Step Manager comparison
349  }
350  }
351 
352 
353  // ----------------------------------------------
354  // Returns the fGhostSafety as the proposedSafety
355  // The SteppingManager will take care of keeping
356  // the smallest one.
357  // ----------------------------------------------
358  return returnedStep;
359 }
ELimited
#define G4ThreadLocal
Definition: tls.hh:89
G4int GetCurrentStepNumber() const
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
G4VPhysicalVolume * GetVolume() const
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
#define DBL_MAX
Definition: templates.hh:83
static void Update(G4FieldTrack *, const G4Track *)

Here is the call graph for this function:

G4VParticleChange * G4FastSimulationManagerProcess::AtRestDoIt ( const G4Track ,
const G4Step  
)
virtual

Implements G4VProcess.

Definition at line 405 of file G4FastSimulationManagerProcess.cc.

406 {
407  return fFastSimulationManager->InvokeAtRestDoIt();
408 }
G4VParticleChange * InvokeAtRestDoIt()

Here is the call graph for this function:

G4double G4FastSimulationManagerProcess::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 378 of file G4FastSimulationManagerProcess.cc.

380 {
381  const G4VPhysicalVolume* currentVolume(0);
382  if ( fIsGhostGeometry ) currentVolume = fPathFinder->GetLocatedVolume(fGhostNavigatorIndex);
383  else currentVolume = track.GetVolume();
384  fFastSimulationManager = currentVolume->GetLogicalVolume()->GetFastSimulationManager();
385  if( fFastSimulationManager )
386  {
387  // Ask for trigger:
388  fFastSimulationTrigger = fFastSimulationManager->AtRestGetFastSimulationManagerTrigger(track, fGhostNavigator);
389  if( fFastSimulationTrigger )
390  {
391  // Dirty trick to take control over stepping. Does anyone will ever use that ?
392  *condition = NotForced;
393  return -1.0;
394  }
395  }
396 
397  // -- no fast simulation occuring there:
398  *condition = NotForced;
399  return DBL_MAX;
400 }
G4double condition(const G4ErrorSymMatrix &m)
G4bool AtRestGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=0)
G4VPhysicalVolume * GetLocatedVolume(G4int navId) const
G4VPhysicalVolume * GetVolume() const
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:

void G4FastSimulationManagerProcess::EndTracking ( )
virtual

Reimplemented from G4VProcess.

Definition at line 229 of file G4FastSimulationManagerProcess.cc.

230 {
231  fIsTrackingTime = false;
232  if ( fIsGhostGeometry ) fTransportationManager->DeActivateNavigator(fGhostNavigator);
233 }
void DeActivateNavigator(G4Navigator *aNavigator)

Here is the call graph for this function:

G4VPhysicalVolume* G4FastSimulationManagerProcess::GetWorldVolume ( ) const
inline

Definition at line 104 of file G4FastSimulationManagerProcess.hh.

104 {return fWorldVolume;}

Here is the caller graph for this function:

G4VParticleChange * G4FastSimulationManagerProcess::PostStepDoIt ( const G4Track ,
const G4Step  
)
virtual

Implements G4VProcess.

Definition at line 280 of file G4FastSimulationManagerProcess.cc.

282 {
283  G4VParticleChange* finalState = fFastSimulationManager->InvokePostStepDoIt();
284 
285  // If the particle is still alive, suspend it to force physics re-initialisation:
286  if (finalState->GetTrackStatus() != fStopAndKill) finalState->ProposeTrackStatus(fSuspend);
287 
288  return finalState;
289 }
G4VParticleChange * InvokePostStepDoIt()
G4TrackStatus GetTrackStatus() const
void ProposeTrackStatus(G4TrackStatus status)

Here is the call graph for this function:

G4double G4FastSimulationManagerProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 241 of file G4FastSimulationManagerProcess.cc.

244 {
245  // -- Get current volume, and check for presence of fast simulation manager.
246  // -- For the case of the navigator for tracking (fGhostNavigatorIndex == 0)
247  // -- we use the track volume. This allows the code to be valid for both
248  // -- cases where the PathFinder is used (G4CoupledTranportation) or not
249  // -- (G4Transportation).
250  const G4VPhysicalVolume* currentVolume(0);
251  if ( fIsGhostGeometry ) currentVolume = fPathFinder->GetLocatedVolume(fGhostNavigatorIndex);
252  else currentVolume = track.GetVolume();
253 
254  if ( currentVolume )
255  {
256  fFastSimulationManager = currentVolume->GetLogicalVolume()->GetFastSimulationManager();
257  if( fFastSimulationManager )
258  {
259  // Ask for trigger:
260  fFastSimulationTrigger = fFastSimulationManager->PostStepGetFastSimulationManagerTrigger(track, fGhostNavigator);
261  if( fFastSimulationTrigger )
262  {
263  // Take control over stepping:
265  return 0.0;
266  }
267  }
268  }
269 
270  // -- no fast simulation occuring there:
271  *condition = NotForced;
272  return DBL_MAX;
273 }
G4double condition(const G4ErrorSymMatrix &m)
G4VPhysicalVolume * GetLocatedVolume(G4int navId) const
G4VPhysicalVolume * GetVolume() const
G4bool PostStepGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=0)
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:

void G4FastSimulationManagerProcess::SetWorldVolume ( G4String  newWorldName)

Definition at line 154 of file G4FastSimulationManagerProcess.cc.

155 {
156  if (fIsTrackingTime)
157  {
159  ed << "G4FastSimulationManagerProcess `" << GetProcessName()
160  << "': changing of world volume at tracking time is not allowed." << G4endl;
161  G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
162  "FastSim002",
163  JustWarning, ed,
164  "Call ignored.");
165  }
166  else
167  {
168  G4VPhysicalVolume* newWorld = fTransportationManager->IsWorldExisting(newWorldName);
169  if (newWorld == 0)
170  {
171  G4ExceptionDescription tellWhatIsWrong;
172  tellWhatIsWrong << "Volume newWorldName = `" << newWorldName
173  << "' is not a parallel world nor the mass world volume."
174  << G4endl;
175  G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
176  "FastSim003",
178  tellWhatIsWrong);
179  }
180  if (verboseLevel>0)
181  {
182  if (fWorldVolume) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
183  << "': changing world volume from '" << fWorldVolume->GetName()
184  << "' to `" << newWorld << "'." << G4endl;
185  else G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
186  << "': setting world volume from to `"<< newWorld->GetName() << "'." << G4endl;
187  }
188  fWorldVolume = newWorld;
189  }
190 }
G4VPhysicalVolume * IsWorldExisting(const G4String &worldName)
G4int verboseLevel
Definition: G4VProcess.hh:368
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

void G4FastSimulationManagerProcess::SetWorldVolume ( G4VPhysicalVolume newWorld)

Definition at line 193 of file G4FastSimulationManagerProcess.cc.

194 {
195  if (newWorld) SetWorldVolume(newWorld->GetName());
196  else
197  {
198  G4ExceptionDescription tellWhatIsWrong;
199  tellWhatIsWrong << "Null pointer passed for world volume." << G4endl;
200  G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4VPhysicalVolume* newWorld)",
201  "FastSim004",
203  tellWhatIsWrong);
204  }
205 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const G4String & GetName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4FastSimulationManagerProcess::StartTracking ( G4Track track)
virtual

Reimplemented from G4VProcess.

Definition at line 211 of file G4FastSimulationManagerProcess.cc.

212 {
213  fIsTrackingTime = true;
214  fIsFirstStep = true;
215 
216  // -- fetch the navigator (and its index) and activate it:
218  fGhostNavigator = transportationManager->GetNavigator(fWorldVolume);
219  fIsGhostGeometry = (fGhostNavigator != transportationManager->GetNavigatorForTracking());
220  if (fIsGhostGeometry) fGhostNavigatorIndex = transportationManager->ActivateNavigator(fGhostNavigator);
221  else fGhostNavigatorIndex = -1;
222 
223  fPathFinder->PrepareNewTrack(track->GetPosition(), track->GetMomentumDirection());
224 }
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=0)
const G4ThreeVector & GetPosition() const
G4Navigator * GetNavigatorForTracking() const
static G4TransportationManager * GetTransportationManager()
G4int ActivateNavigator(G4Navigator *aNavigator)
const G4ThreeVector & GetMomentumDirection() const
G4Navigator * GetNavigator(const G4String &worldName)

Here is the call graph for this function:

void G4FastSimulationManagerProcess::Verbose ( ) const

Definition at line 411 of file G4FastSimulationManagerProcess.cc.

412 {
413  /* G4cout << " >>>>> Trigger Status : ";
414  switch(fFastSimulationManager->GetTriggerStatus())
415  {
416  case NoModel:
417  G4cout << "NoModel" << G4endl;
418  break;
419  case OnBoundaryButLeaving:
420  G4cout << "OnBoundaryButLeaving" << G4endl;
421  break;
422  case OneModelTrigger:
423  G4cout << "OneModelTrigger" << G4endl;
424  break;
425  case NoModelTrigger:
426  G4cout << "NoModelTrigger" << G4endl;
427  break;
428  case Undefined:
429  G4cout << "Undefined" << G4endl;
430  break;
431  default:
432  G4cout << " Bizarre..." << G4endl;
433  break;
434  }*/
435 }

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