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

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.

00055                                                             : 
00056   G4VProcess(processName,theType),
00057   fWorldVolume(0),
00058   fIsTrackingTime(false),
00059   fGhostNavigator(0),
00060   fGhostNavigatorIndex(-1),
00061   fIsGhostGeometry(false),
00062   fGhostSafety(-1.0),
00063   fFieldTrack('0'),
00064   fFastSimulationManager(0),
00065   fFastSimulationTrigger(false)
00066 {
00067   // -- set Process Sub Type
00068   SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
00069 
00070 
00071   fPathFinder            = G4PathFinder::GetInstance();
00072   fTransportationManager = G4TransportationManager::GetTransportationManager();
00073   
00074   SetWorldVolume(fTransportationManager->GetNavigatorForTracking()->GetWorldVolume()->GetName());
00075   if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName() 
00076                              << "' is created, and will message geometry with world volume `" 
00077                              << fWorldVolume->GetName() << "'." << G4endl;
00078   G4GlobalFastSimulationManager::GetGlobalFastSimulationManager()->AddFSMP(this);
00079 }

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.

00085                                                                 :
00086   G4VProcess(processName,theType),
00087   fWorldVolume(0),
00088   fIsTrackingTime(false),
00089   fGhostNavigator(0),
00090   fGhostNavigatorIndex(-1),
00091   fIsGhostGeometry(false),
00092   fGhostSafety(-1.0),
00093   fFieldTrack('0'),
00094   fFastSimulationManager(0),
00095   fFastSimulationTrigger(false)
00096 {
00097   // -- set Process Sub Type
00098   SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
00099 
00100 
00101   fPathFinder            = G4PathFinder::GetInstance();
00102   fTransportationManager = G4TransportationManager::GetTransportationManager();
00103 
00104   SetWorldVolume(worldVolumeName);
00105   if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName() 
00106                              << "' is created, and will message geometry with world volume `" 
00107                              << fWorldVolume->GetName() << "'." << G4endl;
00108   G4GlobalFastSimulationManager::GetGlobalFastSimulationManager()->AddFSMP(this);
00109 }

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.

00115                                                            :
00116   G4VProcess(processName,theType),
00117   fWorldVolume(0),
00118   fIsTrackingTime(false),
00119   fGhostNavigator(0),
00120   fGhostNavigatorIndex(-1),
00121   fIsGhostGeometry(false),
00122   fGhostSafety(-1.0),
00123   fFieldTrack('0'),
00124   fFastSimulationManager(0),
00125   fFastSimulationTrigger(false)
00126 {
00127   // -- set Process Sub Type
00128   SetProcessSubType(static_cast<int>(FASTSIM_ManagerProcess));
00129   
00130 
00131   fPathFinder            = G4PathFinder::GetInstance();
00132   fTransportationManager = G4TransportationManager::GetTransportationManager();
00133 
00134   SetWorldVolume(worldVolume);
00135   if (verboseLevel>0) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName() 
00136                              << "' is created, and will message geometry with world volume `" 
00137                              << fWorldVolume->GetName() << "'." << G4endl;
00138   G4GlobalFastSimulationManager::GetGlobalFastSimulationManager()->AddFSMP(this);
00139 }

G4FastSimulationManagerProcess::~G4FastSimulationManagerProcess (  )  [virtual]

Definition at line 142 of file G4FastSimulationManagerProcess.cc.

References G4GlobalFastSimulationManager::GetGlobalFastSimulationManager(), and G4GlobalFastSimulationManager::RemoveFSMP().


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().

00359 {
00360   fDummyParticleChange.Initialize(track);
00361   return &fDummyParticleChange;
00362 }

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, G4Track::GetCurrentStepNumber(), G4Track::GetVolume(), kDoNot, kSharedOther, kSharedTransport, kUnique, NotCandidateForSelection, and G4FieldTrackUpdator::Update().

00298 {
00299   
00300   *selection            = NotCandidateForSelection;
00301   G4double returnedStep = DBL_MAX;
00302 
00303   // ---------------------------------------------------
00304   // -- Below code valid for ghost geometry, otherwise
00305   // -- useless for fast simulation attached to mass
00306   // -- geometry. Warn user in case along used for 
00307   // -- mass geometry ?
00308   // --------------------------------------------------
00309   if ( fIsGhostGeometry )
00310     {
00311       static G4FieldTrack endTrack('0');
00312       static ELimited eLimited;
00313       
00314       if (previousStepSize > 0.) fGhostSafety -= previousStepSize;
00315       if (fGhostSafety < 0.)     fGhostSafety = 0.0;
00316       
00317       // ------------------------------------------
00318       // Determination of the proposed step length:
00319       // ------------------------------------------
00320       if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
00321         {
00322           // -- No chance to limit the step, as proposed move inside safety
00323           returnedStep   = currentMinimumStep;
00324           proposedSafety = fGhostSafety - currentMinimumStep;
00325         }
00326       else
00327         {
00328           // -- Proposed move exceeds safety, need to state
00329           G4FieldTrackUpdator::Update(&fFieldTrack, &track);
00330           returnedStep = fPathFinder->ComputeStep(fFieldTrack,
00331                                                   currentMinimumStep,
00332                                                   fGhostNavigatorIndex,
00333                                                   track.GetCurrentStepNumber(),
00334                                                   fGhostSafety,
00335                                                   eLimited,
00336                                                   endTrack,
00337                                                   track.GetVolume());
00338           
00339           if(eLimited == kDoNot) fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition()); // -- step no limited by ghost
00340           proposedSafety = fGhostSafety;
00341           if      (eLimited == kUnique || eLimited == kSharedOther) *selection     = CandidateForSelection;
00342           else if (eLimited == kSharedTransport)                     returnedStep *= (1.0 + 1.0e-9);   // -- Expand to disable its selection in Step Manager comparison
00343         }
00344     }
00345   
00346 
00347   // ----------------------------------------------
00348   // Returns the fGhostSafety as the proposedSafety
00349   // The SteppingManager will take care of keeping
00350   // the smallest one.
00351   // ----------------------------------------------
00352   return returnedStep;
00353 }

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

Implements G4VProcess.

Definition at line 399 of file G4FastSimulationManagerProcess.cc.

References G4FastSimulationManager::InvokeAtRestDoIt().

00400 {
00401   return fFastSimulationManager->InvokeAtRestDoIt();
00402 }

G4double G4FastSimulationManagerProcess::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
) [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.

00374 {
00375   const G4VPhysicalVolume* currentVolume(0);
00376   if ( fIsGhostGeometry )  currentVolume = fPathFinder->GetLocatedVolume(fGhostNavigatorIndex);
00377   else                     currentVolume = track.GetVolume();
00378   fFastSimulationManager = currentVolume->GetLogicalVolume()->GetFastSimulationManager();
00379   if( fFastSimulationManager )
00380     {
00381       // Ask for trigger:
00382       fFastSimulationTrigger = fFastSimulationManager->AtRestGetFastSimulationManagerTrigger(track, fGhostNavigator);
00383       if( fFastSimulationTrigger )
00384         {
00385           // Dirty trick to take control over stepping. Does anyone will ever use that ?
00386           *condition = NotForced;
00387           return -1.0;
00388         }
00389     }
00390   
00391   // -- no fast simulation occuring there:
00392   *condition = NotForced;
00393   return DBL_MAX;
00394 }

void G4FastSimulationManagerProcess::EndTracking (  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 228 of file G4FastSimulationManagerProcess.cc.

References G4TransportationManager::DeActivateNavigator().

00229 {
00230   fIsTrackingTime = false;
00231   if ( fIsGhostGeometry ) fTransportationManager->DeActivateNavigator(fGhostNavigator); 
00232 }

G4VPhysicalVolume* G4FastSimulationManagerProcess::GetWorldVolume (  )  const [inline]

Definition at line 104 of file G4FastSimulationManagerProcess.hh.

00104 {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().

00281 {
00282   G4VParticleChange* finalState = fFastSimulationManager->InvokePostStepDoIt();
00283   
00284   // If the particle is still alive, suspend it to force physics re-initialisation:
00285   if (finalState->GetTrackStatus() != fStopAndKill) finalState->ProposeTrackStatus(fSuspend);
00286   
00287   return finalState;
00288 }

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().

00243 {
00244   // -- Get current volume, and check for presence of fast simulation manager.
00245   // -- For the case of the navigator for tracking (fGhostNavigatorIndex == 0)
00246   // -- we use the track volume. This allows the code to be valid for both
00247   // -- cases where the PathFinder is used (G4CoupledTranportation) or not
00248   // -- (G4Transportation).
00249   const G4VPhysicalVolume* currentVolume(0);
00250   if ( fIsGhostGeometry )  currentVolume = fPathFinder->GetLocatedVolume(fGhostNavigatorIndex);
00251   else                     currentVolume = track.GetVolume();
00252 
00253   if ( currentVolume )
00254     {
00255       fFastSimulationManager = currentVolume->GetLogicalVolume()->GetFastSimulationManager();
00256       if( fFastSimulationManager )
00257         {
00258           // Ask for trigger:
00259           fFastSimulationTrigger = fFastSimulationManager->PostStepGetFastSimulationManagerTrigger(track, fGhostNavigator);
00260           if( fFastSimulationTrigger )
00261             {
00262               // Take control over stepping:
00263               *condition = ExclusivelyForced;
00264               return 0.0;
00265             }
00266         }     
00267     }
00268   
00269   // -- no fast simulation occuring there:
00270   *condition = NotForced;
00271   return DBL_MAX;
00272 }

void G4FastSimulationManagerProcess::SetWorldVolume ( G4VPhysicalVolume  ) 

Definition at line 190 of file G4FastSimulationManagerProcess.cc.

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

00191 {
00192   if (newWorld) SetWorldVolume(newWorld->GetName());
00193   else
00194     {
00195       G4ExceptionDescription  tellWhatIsWrong;
00196       tellWhatIsWrong << "Null pointer passed for world volume." << G4endl;
00197       G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4VPhysicalVolume* newWorld)",
00198                   "FastSim004",
00199                   FatalException,
00200                   tellWhatIsWrong); 
00201     }
00202 }

void G4FastSimulationManagerProcess::SetWorldVolume ( G4String   ) 

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().

00152 {
00153   if (fIsTrackingTime)
00154     {
00155       G4ExceptionDescription ed;
00156       ed << "G4FastSimulationManagerProcess `" << GetProcessName()
00157          << "': changing of world volume at tracking time is not allowed." << G4endl;
00158       G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
00159                   "FastSim002",
00160                   JustWarning, ed,
00161                   "Call ignored.");
00162     }
00163   else
00164     {
00165       G4VPhysicalVolume* newWorld = fTransportationManager->IsWorldExisting(newWorldName);
00166       if (newWorld == 0)
00167         {
00168           G4ExceptionDescription  tellWhatIsWrong;
00169           tellWhatIsWrong << "Volume newWorldName = `" <<  newWorldName
00170                           << "' is not a parallel world nor the mass world volume."
00171                           << G4endl;
00172           G4Exception("G4FastSimulationManagerProcess::SetWorldVolume(const G4String)",
00173                       "FastSim003",
00174                       FatalException,
00175                       tellWhatIsWrong);
00176         }
00177       if (verboseLevel>0)
00178       {
00179         if (fWorldVolume) G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
00180                                  << "': changing world volume from '"  << fWorldVolume->GetName() 
00181                                  << "' to `" << newWorld << "'." << G4endl;
00182         else              G4cout << "G4FastSimulationManagerProcess `" << GetProcessName()
00183                                  << "': setting world volume from to `"<< newWorld->GetName() << "'." << G4endl;
00184       }
00185       fWorldVolume = newWorld;
00186     }
00187 }

void G4FastSimulationManagerProcess::StartTracking ( G4Track  )  [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().

00211 {
00212   fIsTrackingTime = true;
00213   fIsFirstStep    = true;
00214   
00215   // -- fetch the navigator (and its index) and activate it:
00216   G4TransportationManager* transportationManager = G4TransportationManager::GetTransportationManager();
00217   fGhostNavigator                            = transportationManager->GetNavigator(fWorldVolume);
00218   fIsGhostGeometry                           = (fGhostNavigator != transportationManager->GetNavigatorForTracking());
00219   if (fIsGhostGeometry) fGhostNavigatorIndex = transportationManager->ActivateNavigator(fGhostNavigator);
00220   else                  fGhostNavigatorIndex = -1;
00221 
00222   fPathFinder->PrepareNewTrack(track->GetPosition(), track->GetMomentumDirection());
00223 }

void G4FastSimulationManagerProcess::Verbose (  )  const

Definition at line 405 of file G4FastSimulationManagerProcess.cc.

00406 {
00407   /*  G4cout << "     >>>>> Trigger Status : ";
00408   switch(fFastSimulationManager->GetTriggerStatus())
00409     {
00410     case NoModel:
00411       G4cout << "NoModel" << G4endl;
00412       break;
00413     case OnBoundaryButLeaving:
00414       G4cout << "OnBoundaryButLeaving" << G4endl;
00415       break;
00416     case OneModelTrigger:
00417       G4cout << "OneModelTrigger" << G4endl;
00418       break;
00419     case NoModelTrigger:
00420       G4cout << "NoModelTrigger" << G4endl;
00421       break;
00422     case Undefined:
00423       G4cout << "Undefined" << G4endl;
00424       break;
00425     default:
00426       G4cout << " Bizarre..." << G4endl;
00427       break;
00428     }*/
00429 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:58 2013 for Geant4 by  doxygen 1.4.7