Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4FastSimulationManagerProcess Class Reference

#include <G4FastSimulationManagerProcess.hh>

Inheritance diagram for G4FastSimulationManagerProcess:
G4VProcess

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.

References G4GlobalFastSimulationManager::AddFSMP(), FASTSIM_ManagerProcess, G4cout, G4endl, G4GlobalFastSimulationManager::GetGlobalFastSimulationManager(), G4PathFinder::GetInstance(), G4VPhysicalVolume::GetName(), G4TransportationManager::GetNavigatorForTracking(), G4VProcess::GetProcessName(), G4TransportationManager::GetTransportationManager(), G4Navigator::GetWorldVolume(), G4VProcess::SetProcessSubType(), SetWorldVolume(), and G4VProcess::verboseLevel.

55  :
56  G4VProcess(processName,theType),
57  fWorldVolume(0),
58  fIsTrackingTime(false),
59  fGhostNavigator(0),
60  fGhostNavigatorIndex(-1),
61  fIsGhostGeometry(false),
62  fGhostSafety(-1.0),
63  fFieldTrack('0'),
64  fFastSimulationManager(0),
65  fFastSimulationTrigger(false)
66 {
67  // -- set Process Sub Type
68  SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
69 
70 
71  fPathFinder = G4PathFinder::GetInstance();
72  fTransportationManager = G4TransportationManager::GetTransportationManager();
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 }
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 *)
G4FastSimulationManagerProcess::G4FastSimulationManagerProcess ( const G4String processName,
const G4String worldVolumeName,
G4ProcessType  theType = fParameterisation 
)

Definition at line 83 of file G4FastSimulationManagerProcess.cc.

References G4GlobalFastSimulationManager::AddFSMP(), FASTSIM_ManagerProcess, G4cout, G4endl, G4GlobalFastSimulationManager::GetGlobalFastSimulationManager(), G4PathFinder::GetInstance(), G4VPhysicalVolume::GetName(), G4VProcess::GetProcessName(), G4TransportationManager::GetTransportationManager(), G4VProcess::SetProcessSubType(), SetWorldVolume(), and G4VProcess::verboseLevel.

85  :
86  G4VProcess(processName,theType),
87  fWorldVolume(0),
88  fIsTrackingTime(false),
89  fGhostNavigator(0),
90  fGhostNavigatorIndex(-1),
91  fIsGhostGeometry(false),
92  fGhostSafety(-1.0),
93  fFieldTrack('0'),
94  fFastSimulationManager(0),
95  fFastSimulationTrigger(false)
96 {
97  // -- set Process Sub Type
98  SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
99 
100 
101  fPathFinder = G4PathFinder::GetInstance();
102  fTransportationManager = G4TransportationManager::GetTransportationManager();
103 
104  SetWorldVolume(worldVolumeName);
105  if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
106  << "' is created, and will message geometry with world volume `"
107  << fWorldVolume->GetName() << "'." << G4endl;
109 }
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 *)
G4FastSimulationManagerProcess::G4FastSimulationManagerProcess ( const G4String processName,
G4VPhysicalVolume worldVolume,
G4ProcessType  theType = fParameterisation 
)

Definition at line 113 of file G4FastSimulationManagerProcess.cc.

References G4GlobalFastSimulationManager::AddFSMP(), FASTSIM_ManagerProcess, G4cout, G4endl, G4GlobalFastSimulationManager::GetGlobalFastSimulationManager(), G4PathFinder::GetInstance(), G4VPhysicalVolume::GetName(), G4VProcess::GetProcessName(), G4TransportationManager::GetTransportationManager(), G4VProcess::SetProcessSubType(), SetWorldVolume(), and G4VProcess::verboseLevel.

115  :
116  G4VProcess(processName,theType),
117  fWorldVolume(0),
118  fIsTrackingTime(false),
119  fGhostNavigator(0),
120  fGhostNavigatorIndex(-1),
121  fIsGhostGeometry(false),
122  fGhostSafety(-1.0),
123  fFieldTrack('0'),
124  fFastSimulationManager(0),
125  fFastSimulationTrigger(false)
126 {
127  // -- set Process Sub Type
128  SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
129 
130 
131  fPathFinder = G4PathFinder::GetInstance();
132  fTransportationManager = G4TransportationManager::GetTransportationManager();
133 
134  SetWorldVolume(worldVolume);
135  if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
136  << "' is created, and will message geometry with world volume `"
137  << fWorldVolume->GetName() << "'." << G4endl;
139 }
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 *)
G4FastSimulationManagerProcess::~G4FastSimulationManagerProcess ( )
virtual

Member Function Documentation

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

Implements G4VProcess.

Definition at line 357 of file G4FastSimulationManagerProcess.cc.

References G4VParticleChange::Initialize().

359 {
360  fDummyParticleChange.Initialize(track);
361  return &fDummyParticleChange;
362 }
virtual void Initialize(const G4Track &)
G4double G4FastSimulationManagerProcess::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
)
virtual

Implements G4VProcess.

Definition at line 293 of file G4FastSimulationManagerProcess.cc.

References CandidateForSelection, G4Navigator::ComputeSafety(), G4PathFinder::ComputeStep(), DBL_MAX, G4ThreadLocal, G4Track::GetCurrentStepNumber(), G4Track::GetVolume(), kDoNot, kSharedOther, kSharedTransport, kUnique, NotCandidateForSelection, and G4FieldTrackUpdator::Update().

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 ; if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ; G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
312  static G4ThreadLocal ELimited *eLimited_G4MT_TLS_ = 0 ; if (!eLimited_G4MT_TLS_) eLimited_G4MT_TLS_ = new ELimited ; ELimited &eLimited = *eLimited_G4MT_TLS_;
313 
314  if (previousStepSize > 0.) fGhostSafety -= previousStepSize;
315  if (fGhostSafety < 0.) fGhostSafety = 0.0;
316 
317  // ------------------------------------------
318  // Determination of the proposed step length:
319  // ------------------------------------------
320  if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
321  {
322  // -- No chance to limit the step, as proposed move inside safety
323  returnedStep = currentMinimumStep;
324  proposedSafety = fGhostSafety - currentMinimumStep;
325  }
326  else
327  {
328  // -- Proposed move exceeds safety, need to state
329  G4FieldTrackUpdator::Update(&fFieldTrack, &track);
330  returnedStep = fPathFinder->ComputeStep(fFieldTrack,
331  currentMinimumStep,
332  fGhostNavigatorIndex,
333  track.GetCurrentStepNumber(),
334  fGhostSafety,
335  eLimited,
336  endTrack,
337  track.GetVolume());
338 
339  if(eLimited == kDoNot) fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition()); // -- step no limited by ghost
340  proposedSafety = fGhostSafety;
341  if (eLimited == kUnique || eLimited == kSharedOther) *selection = CandidateForSelection;
342  else if (eLimited == kSharedTransport) returnedStep *= (1.0 + 1.0e-9); // -- Expand to disable its selection in Step Manager comparison
343  }
344  }
345 
346 
347  // ----------------------------------------------
348  // Returns the fGhostSafety as the proposedSafety
349  // The SteppingManager will take care of keeping
350  // the smallest one.
351  // ----------------------------------------------
352  return returnedStep;
353 }
ELimited
#define G4ThreadLocal
Definition: tls.hh:52
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 *)
G4VParticleChange * G4FastSimulationManagerProcess::AtRestDoIt ( const G4Track ,
const G4Step  
)
virtual

Implements G4VProcess.

Definition at line 399 of file G4FastSimulationManagerProcess.cc.

References G4FastSimulationManager::InvokeAtRestDoIt().

400 {
401  return fFastSimulationManager->InvokeAtRestDoIt();
402 }
G4VParticleChange * InvokeAtRestDoIt()
G4double G4FastSimulationManagerProcess::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 372 of file G4FastSimulationManagerProcess.cc.

References G4FastSimulationManager::AtRestGetFastSimulationManagerTrigger(), DBL_MAX, G4LogicalVolume::GetFastSimulationManager(), G4PathFinder::GetLocatedVolume(), G4VPhysicalVolume::GetLogicalVolume(), G4Track::GetVolume(), and NotForced.

374 {
375  const G4VPhysicalVolume* currentVolume(0);
376  if ( fIsGhostGeometry ) currentVolume = fPathFinder->GetLocatedVolume(fGhostNavigatorIndex);
377  else currentVolume = track.GetVolume();
378  fFastSimulationManager = currentVolume->GetLogicalVolume()->GetFastSimulationManager();
379  if( fFastSimulationManager )
380  {
381  // Ask for trigger:
382  fFastSimulationTrigger = fFastSimulationManager->AtRestGetFastSimulationManagerTrigger(track, fGhostNavigator);
383  if( fFastSimulationTrigger )
384  {
385  // Dirty trick to take control over stepping. Does anyone will ever use that ?
386  *condition = NotForced;
387  return -1.0;
388  }
389  }
390 
391  // -- no fast simulation occuring there:
392  *condition = NotForced;
393  return DBL_MAX;
394 }
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
void G4FastSimulationManagerProcess::EndTracking ( )
virtual

Reimplemented from G4VProcess.

Definition at line 228 of file G4FastSimulationManagerProcess.cc.

References G4TransportationManager::DeActivateNavigator().

229 {
230  fIsTrackingTime = false;
231  if ( fIsGhostGeometry ) fTransportationManager->DeActivateNavigator(fGhostNavigator);
232 }
void DeActivateNavigator(G4Navigator *aNavigator)
G4VPhysicalVolume* G4FastSimulationManagerProcess::GetWorldVolume ( ) const
inline

Definition at line 104 of file G4FastSimulationManagerProcess.hh.

104 {return fWorldVolume;}
G4VParticleChange * G4FastSimulationManagerProcess::PostStepDoIt ( const G4Track ,
const G4Step  
)
virtual

Implements G4VProcess.

Definition at line 279 of file G4FastSimulationManagerProcess.cc.

References fStopAndKill, fSuspend, G4VParticleChange::GetTrackStatus(), G4FastSimulationManager::InvokePostStepDoIt(), and G4VParticleChange::ProposeTrackStatus().

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 }
G4VParticleChange * InvokePostStepDoIt()
G4TrackStatus GetTrackStatus() const
void ProposeTrackStatus(G4TrackStatus status)
G4double G4FastSimulationManagerProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 240 of file G4FastSimulationManagerProcess.cc.

References DBL_MAX, ExclusivelyForced, G4LogicalVolume::GetFastSimulationManager(), G4PathFinder::GetLocatedVolume(), G4VPhysicalVolume::GetLogicalVolume(), G4Track::GetVolume(), NotForced, and G4FastSimulationManager::PostStepGetFastSimulationManagerTrigger().

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:
270  *condition = NotForced;
271  return DBL_MAX;
272 }
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
void G4FastSimulationManagerProcess::SetWorldVolume ( G4String  newWorldName)

Definition at line 151 of file G4FastSimulationManagerProcess.cc.

References FatalException, G4cout, G4endl, G4Exception(), G4VPhysicalVolume::GetName(), G4VProcess::GetProcessName(), G4TransportationManager::IsWorldExisting(), JustWarning, and G4VProcess::verboseLevel.

Referenced by G4FastSimulationManagerProcess(), and SetWorldVolume().

152 {
153  if (fIsTrackingTime)
154  {
156  ed << "G4FastSimulationManagerProcess `" << GetProcessName()
157  << "': changing of world volume at tracking time is not allowed." << G4endl;
158  G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
159  "FastSim002",
160  JustWarning, ed,
161  "Call ignored.");
162  }
163  else
164  {
165  G4VPhysicalVolume* newWorld = fTransportationManager->IsWorldExisting(newWorldName);
166  if (newWorld == 0)
167  {
168  G4ExceptionDescription tellWhatIsWrong;
169  tellWhatIsWrong << "Volume newWorldName = `" << newWorldName
170  << "' is not a parallel world nor the mass world volume."
171  << G4endl;
172  G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
173  "FastSim003",
175  tellWhatIsWrong);
176  }
177  if (verboseLevel>0)
178  {
179  if (fWorldVolume) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
180  << "': changing world volume from '" << fWorldVolume->GetName()
181  << "' to `" << newWorld << "'." << G4endl;
182  else G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
183  << "': setting world volume from to `"<< newWorld->GetName() << "'." << G4endl;
184  }
185  fWorldVolume = newWorld;
186  }
187 }
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
void G4FastSimulationManagerProcess::SetWorldVolume ( G4VPhysicalVolume newWorld)

Definition at line 190 of file G4FastSimulationManagerProcess.cc.

References FatalException, G4endl, G4Exception(), G4VPhysicalVolume::GetName(), and SetWorldVolume().

191 {
192  if (newWorld) SetWorldVolume(newWorld->GetName());
193  else
194  {
195  G4ExceptionDescription tellWhatIsWrong;
196  tellWhatIsWrong << "Null pointer passed for world volume." << G4endl;
197  G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4VPhysicalVolume* newWorld)",
198  "FastSim004",
200  tellWhatIsWrong);
201  }
202 }
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
void G4FastSimulationManagerProcess::StartTracking ( G4Track track)
virtual

Reimplemented from G4VProcess.

Definition at line 210 of file G4FastSimulationManagerProcess.cc.

References G4TransportationManager::ActivateNavigator(), G4Track::GetMomentumDirection(), G4TransportationManager::GetNavigator(), G4TransportationManager::GetNavigatorForTracking(), G4Track::GetPosition(), G4TransportationManager::GetTransportationManager(), and G4PathFinder::PrepareNewTrack().

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 }
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)
void G4FastSimulationManagerProcess::Verbose ( ) const

Definition at line 405 of file G4FastSimulationManagerProcess.cc.

406 {
407  /* G4cout << " >>>>> Trigger Status : ";
408  switch(fFastSimulationManager->GetTriggerStatus())
409  {
410  case NoModel:
411  G4cout << "NoModel" << G4endl;
412  break;
413  case OnBoundaryButLeaving:
414  G4cout << "OnBoundaryButLeaving" << G4endl;
415  break;
416  case OneModelTrigger:
417  G4cout << "OneModelTrigger" << G4endl;
418  break;
419  case NoModelTrigger:
420  G4cout << "NoModelTrigger" << G4endl;
421  break;
422  case Undefined:
423  G4cout << "Undefined" << G4endl;
424  break;
425  default:
426  G4cout << " Bizarre..." << G4endl;
427  break;
428  }*/
429 }

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