Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FastSimulationManagerProcess Class Reference

#include <G4FastSimulationManagerProcess.hh>

+ Inheritance 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)
 
 ~G4FastSimulationManagerProcess () override
 
G4VPhysicalVolumeGetWorldVolume () const
 
void SetWorldVolume (G4String)
 
void SetWorldVolume (G4VPhysicalVolume *)
 
void StartTracking (G4Track *) override
 
void EndTracking () override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection) override
 
G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &step) override
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *) override
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &) override
 
void Verbose () const
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool 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 const G4VProcessGetCreatorProcess () const
 
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
 
virtual void ProcessDescription (std::ostream &outfile) 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 prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 75 of file G4FastSimulationManagerProcess.hh.

Constructor & Destructor Documentation

◆ G4FastSimulationManagerProcess() [1/3]

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

Definition at line 51 of file G4FastSimulationManagerProcess.cc.

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}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void AddFSMP(G4FastSimulationManagerProcess *)
static G4GlobalFastSimulationManager * GetGlobalFastSimulationManager()
G4VPhysicalVolume * GetWorldVolume() const
static G4PathFinder * GetInstance()
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
const G4String & GetName() const
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const

◆ G4FastSimulationManagerProcess() [2/3]

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

Definition at line 79 of file G4FastSimulationManagerProcess.cc.

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}

◆ G4FastSimulationManagerProcess() [3/3]

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

Definition at line 108 of file G4FastSimulationManagerProcess.cc.

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}

◆ ~G4FastSimulationManagerProcess()

G4FastSimulationManagerProcess::~G4FastSimulationManagerProcess ( )
override

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Definition at line 325 of file G4FastSimulationManagerProcess.cc.

327{
328 fDummyParticleChange.Initialize(track);
329 return &fDummyParticleChange;
330}
virtual void Initialize(const G4Track &)

◆ AlongStepGetPhysicalInteractionLength()

G4double G4FastSimulationManagerProcess::AlongStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & proposedSafety,
G4GPILSelection * selection )
overridevirtual

Implements G4VProcess.

Definition at line 265 of file G4FastSimulationManagerProcess.cc.

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}
@ CandidateForSelection
@ NotCandidateForSelection
@ kDoNot
@ kUnique
@ kSharedOther
@ kSharedTransport
double G4double
Definition G4Types.hh:83
static void Update(G4FieldTrack *, const G4Track *)
G4ThreeVector GetPosition() const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
G4VPhysicalVolume * GetVolume() const
G4int GetCurrentStepNumber() const
#define DBL_MAX
Definition templates.hh:62
#define G4ThreadLocal
Definition tls.hh:77

◆ AtRestDoIt()

G4VParticleChange * G4FastSimulationManagerProcess::AtRestDoIt ( const G4Track & ,
const G4Step &  )
overridevirtual

Implements G4VProcess.

Definition at line 366 of file G4FastSimulationManagerProcess.cc.

367{
368 return fFastSimulationManager->InvokeAtRestDoIt();
369}
G4VParticleChange * InvokeAtRestDoIt()

◆ AtRestGetPhysicalInteractionLength()

G4double G4FastSimulationManagerProcess::AtRestGetPhysicalInteractionLength ( const G4Track & track,
G4ForceCondition * condition )
overridevirtual

Implements G4VProcess.

Definition at line 338 of file G4FastSimulationManagerProcess.cc.

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}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
G4bool AtRestGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=nullptr)
G4VPhysicalVolume * GetLocatedVolume(G4int navId) const

◆ EndTracking()

void G4FastSimulationManagerProcess::EndTracking ( )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 209 of file G4FastSimulationManagerProcess.cc.

210{
211 fIsTrackingTime = false;
212 if (fIsGhostGeometry) fTransportationManager->DeActivateNavigator(fGhostNavigator);
213}
void DeActivateNavigator(G4Navigator *aNavigator)

◆ GetWorldVolume()

G4VPhysicalVolume * G4FastSimulationManagerProcess::GetWorldVolume ( ) const
inline

Definition at line 95 of file G4FastSimulationManagerProcess.hh.

95{ return fWorldVolume; }

Referenced by G4FastSimulationPhysics::ConstructProcess().

◆ PostStepDoIt()

G4VParticleChange * G4FastSimulationManagerProcess::PostStepDoIt ( const G4Track & ,
const G4Step &  )
overridevirtual

Implements G4VProcess.

Definition at line 255 of file G4FastSimulationManagerProcess.cc.

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}
@ fSuspend
@ fStopAndKill
G4VParticleChange * InvokePostStepDoIt()
void ProposeTrackStatus(G4TrackStatus status)
G4TrackStatus GetTrackStatus() const

◆ PostStepGetPhysicalInteractionLength()

G4double G4FastSimulationManagerProcess::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Implements G4VProcess.

Definition at line 219 of file G4FastSimulationManagerProcess.cc.

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}
@ ExclusivelyForced
G4bool PostStepGetFastSimulationManagerTrigger(const G4Track &, const G4Navigator *a=nullptr)

◆ SetWorldVolume() [1/2]

void G4FastSimulationManagerProcess::SetWorldVolume ( G4String newWorldName)

Definition at line 145 of file G4FastSimulationManagerProcess.cc.

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}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4VPhysicalVolume * IsWorldExisting(const G4String &worldName)

Referenced by G4FastSimulationManagerProcess(), G4FastSimulationManagerProcess(), G4FastSimulationManagerProcess(), and SetWorldVolume().

◆ SetWorldVolume() [2/2]

void G4FastSimulationManagerProcess::SetWorldVolume ( G4VPhysicalVolume * newWorld)

Definition at line 176 of file G4FastSimulationManagerProcess.cc.

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}

◆ StartTracking()

void G4FastSimulationManagerProcess::StartTracking ( G4Track * track)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 191 of file G4FastSimulationManagerProcess.cc.

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}
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=nullptr)
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
G4int ActivateNavigator(G4Navigator *aNavigator)
G4Navigator * GetNavigator(const G4String &worldName)

◆ Verbose()

void G4FastSimulationManagerProcess::Verbose ( ) const
inline

Definition at line 128 of file G4FastSimulationManagerProcess.hh.

128: will be remove in next major release")]] void Verbose() const {}

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