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

#include <G4ParallelGeometriesLimiterProcess.hh>

+ Inheritance diagram for G4ParallelGeometriesLimiterProcess:

Public Member Functions

 G4ParallelGeometriesLimiterProcess (const G4String &processName="biasLimiter")
 
virtual ~G4ParallelGeometriesLimiterProcess ()=default
 
void AddParallelWorld (const G4String &parallelWorldName)
 
void RemoveParallelWorld (const G4String &parallelWorldName)
 
const std::vector< G4VPhysicalVolume * > & GetParallelWorlds () const
 
G4int GetParallelWorldIndex (const G4VPhysicalVolume *parallelWorld) const
 
G4int GetParallelWorldIndex (const G4String &parallelWorldName) const
 
const std::vector< G4Navigator * > & GetActiveNavigators () const
 
const G4NavigatorGetNavigator (G4int worldIndex) const
 
const std::vector< const G4VPhysicalVolume * > & GetCurrentVolumes () const
 
const std::vector< const G4VPhysicalVolume * > & GetPreviousVolumes () const
 
const G4VPhysicalVolumeGetCurrentVolume (G4int worldIndex) const
 
const G4VPhysicalVolumeGetPreviousVolume (G4int worldIndex) const
 
const std::vector< G4bool > & GetIsLimiting () const
 
const std::vector< G4bool > & GetWasLimiting () const
 
G4bool GetIsLimiting (G4int worldIndex) const
 
G4bool GetWasLimiting (G4int worldIndex) const
 
void StartTracking (G4Track *)
 
void EndTracking ()
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *)
 
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 &)
 
virtual void SetProcessManager (const G4ProcessManager *)
 
- 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 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 49 of file G4ParallelGeometriesLimiterProcess.hh.

Constructor & Destructor Documentation

◆ G4ParallelGeometriesLimiterProcess()

G4ParallelGeometriesLimiterProcess::G4ParallelGeometriesLimiterProcess ( const G4String & processName = "biasLimiter")

Definition at line 39 of file G4ParallelGeometriesLimiterProcess.cc.

41 : G4VProcess(processName, fParallel)
42{
43 // -- Process Sub Type ?
44
45 fPathFinder = G4PathFinder::GetInstance();
47}
@ fParallel
static G4PathFinder * GetInstance()
static G4TransportationManager * GetTransportationManager()
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition G4VProcess.cc:46

◆ ~G4ParallelGeometriesLimiterProcess()

virtual G4ParallelGeometriesLimiterProcess::~G4ParallelGeometriesLimiterProcess ( )
virtualdefault

Member Function Documentation

◆ AddParallelWorld()

void G4ParallelGeometriesLimiterProcess::AddParallelWorld ( const G4String & parallelWorldName)

Definition at line 52 of file G4ParallelGeometriesLimiterProcess.cc.

54{
55 // -- Refuse adding parallel geometry during tracking time:
56 if (fIsTrackingTime)
57 {
59 ed << "G4ParallelGeometriesLimiterProcess `" << GetProcessName()
60 << "': adding a parallel world volume at tracking time is not allowed."
61 << G4endl;
62 G4Exception("G4ParallelGeometriesLimiterProcess::AddParallelWorld(..)",
63 "BIAS.GEN.21", JustWarning, ed, "Call ignored.");
64 return;
65 }
66 else
67 {
68 G4VPhysicalVolume* newWorld = fTransportationManager->IsWorldExisting( parallelWorldName );
69
70 // -- Fatal exception if requested world does not exist:
71 if (newWorld == nullptr)
72 {
73 G4ExceptionDescription tellWhatIsWrong;
74 tellWhatIsWrong << "Volume `" << parallelWorldName
75 << "' is not a parallel world nor the mass world volume."
76 << G4endl;
77 G4Exception("G4ParallelGeometriesLimiterProcess::SetWorldVolume(..)",
78 "BIAS.GEN.22", FatalException, tellWhatIsWrong);
79 }
80
81 // -- Protection against adding the mass geometry world as parallel world:
82 if ( newWorld == fTransportationManager->GetNavigatorForTracking()->GetWorldVolume() )
83 {
85 ed << "G4ParallelGeometriesLimiterProcess `" << GetProcessName()
86 << "': trying to add the world volume for tracking as a parallel world."
87 << G4endl;
88 G4Exception("G4ParallelGeometriesLimiterProcess::AddParallelWorld(..)",
89 "BIAS.GEN.23", JustWarning, ed, "Call ignored.");
90 return;
91 }
92
93 // -- Add parallel world, taking care it is not in the list yet:
94 G4bool isNew = true;
95 for ( auto knownWorld : fParallelWorlds )
96 {
97 if ( knownWorld == newWorld ) { isNew = false; }
98 }
99 if ( isNew )
100 {
101 fParallelWorlds.push_back( newWorld );
102 }
103 else
104 {
106 ed << "G4ParallelGeometriesLimiterProcess `" << GetProcessName()
107 << "': trying to re-add the parallel world volume `"
108 << parallelWorldName << "'." << G4endl;
109 G4Exception("G4ParallelGeometriesLimiterProcess::AddParallelWorld(..)",
110 "BIAS.GEN.24", JustWarning, ed, "Call ignored.");
111 return;
112 }
113 }
114}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
const G4String & GetProcessName() const

◆ AlongStepDoIt()

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

Implements G4VProcess.

Definition at line 354 of file G4ParallelGeometriesLimiterProcess.cc.

356{
357 fDummyParticleChange.Initialize(track);
358 return &fDummyParticleChange;
359}

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 236 of file G4ParallelGeometriesLimiterProcess.cc.

242{
243 // -- Init:
244 // -- Note that the returnedStep must be physically meaningful,
245 // -- even if we return NotCandidateForSelection as condition;
246 // -- the reason is that the stepping manager always takes the smallest
247 // -- alongstep among the returned ones (point related
248 // -- to geometry step length wrt to true path length).
249 *selection = NotCandidateForSelection;
250 G4double returnedStep = DBL_MAX;
251
252 // -- G4FieldTrack and ELimited:
253 static G4ThreadLocal G4FieldTrack* endTrack_MT = nullptr;
254 if (!endTrack_MT) endTrack_MT = new G4FieldTrack ('0');
255 G4FieldTrack& endTrack = *endTrack_MT;
256
257 static G4ThreadLocal ELimited* eLimited_MT = nullptr;
258 if (!eLimited_MT) eLimited_MT = new ELimited;
259 ELimited &eLimited = *eLimited_MT;
260
261 // -------------------
262 // -- Update safeties:
263 // -------------------
264 if ( previousStepSize > 0.0 )
265 {
266 for ( auto& parallelWorldSafety : fParallelWorldSafeties )
267 {
268 parallelWorldSafety -= previousStepSize;
269 if ( parallelWorldSafety < 0. ) { parallelWorldSafety = 0.0; }
270 fParallelWorldSafety = parallelWorldSafety < fParallelWorldSafety
271 ? parallelWorldSafety : fParallelWorldSafety;
272 }
273 }
274
275 // ------------------------------------------
276 // Determination of the proposed step length:
277 // ------------------------------------------
278 if ( ( currentMinimumStep <= fParallelWorldSafety )
279 && ( currentMinimumStep > 0. ) )
280 {
281 // -- No chance to limit the step, as proposed move inside safety
282
283 returnedStep = currentMinimumStep;
284 proposedSafety = fParallelWorldSafety - currentMinimumStep;
285 }
286 else
287 {
288 // -- Proposed move exceeds common safety, need to state
289 G4double smallestReturnedStep = -1.0;
290 ELimited eLimitedForSmallestStep = kDoNot;
291 for ( std::size_t i = 0 ; i < fParallelWorldNavigatorIndeces.size() ; ++i )
292 {
293 // -- Update safety of geometries having safety smaller than current minimum step
294 if ( currentMinimumStep >= fParallelWorldSafeties[i] )
295 {
296 G4FieldTrackUpdator::Update(&fFieldTrack, &track);
297 G4double tmpReturnedStep = fPathFinder->ComputeStep(fFieldTrack,
298 currentMinimumStep,
299 fParallelWorldNavigatorIndeces[i],
300 track.GetCurrentStepNumber(),
301 fParallelWorldSafeties[i],
302 eLimited, endTrack, track.GetVolume());
303
304 if ( ( smallestReturnedStep < 0.0 ) || ( tmpReturnedStep <= smallestReturnedStep ) )
305 {
306 smallestReturnedStep = tmpReturnedStep;
307 eLimitedForSmallestStep = eLimited;
308 }
309
310 if (eLimited == kDoNot)
311 {
312 // -- Step not limited by this geometry
313 fParallelWorldSafeties[i] = fParallelWorldNavigators[i]->ComputeSafety(endTrack.GetPosition());
314 fParallelWorldIsLimiting[i] = false;
315 }
316 else
317 {
318 fParallelWorldIsLimiting[i] = true;
319 }
320 }
321
322 // -- update with smallest safety:
323 fParallelWorldSafety = fParallelWorldSafeties[i] < fParallelWorldSafety
324 ? fParallelWorldSafeties[i] : fParallelWorldSafety;
325 }
326
327 // -- no geometry limitation among all geometries, can return currentMinimumStep (or DBL_MAX):
328 // -- Beware : the returnedStep must be physically meaningful, even if we say "NotCandidateForSelection" !
329 if ( eLimitedForSmallestStep == kDoNot )
330 {
331 returnedStep = currentMinimumStep;
332 }
333 // -- proposed step length of limiting geometry:
334 if ( eLimitedForSmallestStep == kUnique ||
335 eLimitedForSmallestStep == kSharedOther )
336 {
337 *selection = CandidateForSelection;
338 returnedStep = smallestReturnedStep;
339 }
340 else if ( eLimitedForSmallestStep == kSharedTransport )
341 {
342 // -- Expand to disable its selection in Step Manager comparison
343 returnedStep = smallestReturnedStep* (1.0 + 1.0e-9);
344 }
345
346 // -- and smallest safety among geometries:
347 proposedSafety = fParallelWorldSafety ;
348 }
349
350 // -- returns step length, and proposedSafety
351 return returnedStep;
352}
@ CandidateForSelection
@ NotCandidateForSelection
@ kDoNot
@ kUnique
@ kSharedOther
@ kSharedTransport
double G4double
Definition G4Types.hh:83
static void Update(G4FieldTrack *, const G4Track *)
G4ThreeVector GetPosition() const
G4VPhysicalVolume * GetVolume() const
G4int GetCurrentStepNumber() const
#define DBL_MAX
Definition templates.hh:62
#define G4ThreadLocal
Definition tls.hh:77

◆ AtRestDoIt()

G4VParticleChange * G4ParallelGeometriesLimiterProcess::AtRestDoIt ( const G4Track & ,
const G4Step &  )
inlinevirtual

Implements G4VProcess.

Definition at line 147 of file G4ParallelGeometriesLimiterProcess.hh.

148 { return nullptr; }

◆ AtRestGetPhysicalInteractionLength()

G4double G4ParallelGeometriesLimiterProcess::AtRestGetPhysicalInteractionLength ( const G4Track & ,
G4ForceCondition *  )
inlinevirtual

Implements G4VProcess.

Definition at line 145 of file G4ParallelGeometriesLimiterProcess.hh.

146 { return DBL_MAX; }

◆ EndTracking()

void G4ParallelGeometriesLimiterProcess::EndTracking ( )
virtual

Reimplemented from G4VProcess.

Definition at line 207 of file G4ParallelGeometriesLimiterProcess.cc.

208{
209 fIsTrackingTime = false;
210 for ( auto parallelWorldNavigator : fParallelWorldNavigators )
211 {
212 fTransportationManager->DeActivateNavigator( parallelWorldNavigator );
213 }
214}

◆ GetActiveNavigators()

const std::vector< G4Navigator * > & G4ParallelGeometriesLimiterProcess::GetActiveNavigators ( ) const
inline

Definition at line 81 of file G4ParallelGeometriesLimiterProcess.hh.

82 { return fParallelWorldNavigators; }

◆ GetCurrentVolume()

const G4VPhysicalVolume * G4ParallelGeometriesLimiterProcess::GetCurrentVolume ( G4int worldIndex) const
inline

Definition at line 103 of file G4ParallelGeometriesLimiterProcess.hh.

104 { return fCurrentVolumes[std::size_t(worldIndex)]; }

◆ GetCurrentVolumes()

const std::vector< const G4VPhysicalVolume * > & G4ParallelGeometriesLimiterProcess::GetCurrentVolumes ( ) const
inline

Definition at line 97 of file G4ParallelGeometriesLimiterProcess.hh.

98 { return fCurrentVolumes; }

◆ GetIsLimiting() [1/2]

const std::vector< G4bool > & G4ParallelGeometriesLimiterProcess::GetIsLimiting ( ) const
inline

Definition at line 108 of file G4ParallelGeometriesLimiterProcess.hh.

109 { return fParallelWorldIsLimiting; }

◆ GetIsLimiting() [2/2]

G4bool G4ParallelGeometriesLimiterProcess::GetIsLimiting ( G4int worldIndex) const
inline

Definition at line 114 of file G4ParallelGeometriesLimiterProcess.hh.

115 { return fParallelWorldIsLimiting[std::size_t(worldIndex)]; }

◆ GetNavigator()

const G4Navigator * G4ParallelGeometriesLimiterProcess::GetNavigator ( G4int worldIndex) const
inline

Definition at line 85 of file G4ParallelGeometriesLimiterProcess.hh.

86 { return fParallelWorldNavigators[std::size_t(worldIndex)]; }

◆ GetParallelWorldIndex() [1/2]

G4int G4ParallelGeometriesLimiterProcess::GetParallelWorldIndex ( const G4String & parallelWorldName) const

Definition at line 411 of file G4ParallelGeometriesLimiterProcess.cc.

413{
414 G4VPhysicalVolume* aWorld = fTransportationManager->IsWorldExisting( parallelWorldName );
415 // note aWorld might be nullptr
416 return GetParallelWorldIndex( aWorld );
417}
G4int GetParallelWorldIndex(const G4VPhysicalVolume *parallelWorld) const

◆ GetParallelWorldIndex() [2/2]

G4int G4ParallelGeometriesLimiterProcess::GetParallelWorldIndex ( const G4VPhysicalVolume * parallelWorld) const

Definition at line 394 of file G4ParallelGeometriesLimiterProcess.cc.

396{
397 G4int toReturn = -1;
398 G4int iWorld = 0;
399 for ( auto world : fParallelWorlds )
400 {
401 if ( world == parallelWorld )
402 {
403 toReturn = iWorld;
404 break;
405 }
406 ++iWorld;
407 }
408 return toReturn;
409}
int G4int
Definition G4Types.hh:85

Referenced by GetParallelWorldIndex().

◆ GetParallelWorlds()

const std::vector< G4VPhysicalVolume * > & G4ParallelGeometriesLimiterProcess::GetParallelWorlds ( ) const
inline

Definition at line 67 of file G4ParallelGeometriesLimiterProcess.hh.

68 { return fParallelWorlds; }

◆ GetPreviousVolume()

const G4VPhysicalVolume * G4ParallelGeometriesLimiterProcess::GetPreviousVolume ( G4int worldIndex) const
inline

Definition at line 105 of file G4ParallelGeometriesLimiterProcess.hh.

106 { return fPreviousVolumes[std::size_t(worldIndex)]; }

◆ GetPreviousVolumes()

const std::vector< const G4VPhysicalVolume * > & G4ParallelGeometriesLimiterProcess::GetPreviousVolumes ( ) const
inline

Definition at line 99 of file G4ParallelGeometriesLimiterProcess.hh.

100 { return fPreviousVolumes; }

◆ GetWasLimiting() [1/2]

const std::vector< G4bool > & G4ParallelGeometriesLimiterProcess::GetWasLimiting ( ) const
inline

Definition at line 110 of file G4ParallelGeometriesLimiterProcess.hh.

111 { return fParallelWorldWasLimiting; }

◆ GetWasLimiting() [2/2]

G4bool G4ParallelGeometriesLimiterProcess::GetWasLimiting ( G4int worldIndex) const
inline

Definition at line 116 of file G4ParallelGeometriesLimiterProcess.hh.

117 { return fParallelWorldWasLimiting[std::size_t(worldIndex)]; }

◆ PostStepDoIt()

G4VParticleChange * G4ParallelGeometriesLimiterProcess::PostStepDoIt ( const G4Track & ,
const G4Step &  )
inlinevirtual

Implements G4VProcess.

Definition at line 132 of file G4ParallelGeometriesLimiterProcess.hh.

132{ return nullptr; }

◆ PostStepGetPhysicalInteractionLength()

G4double G4ParallelGeometriesLimiterProcess::PostStepGetPhysicalInteractionLength ( const G4Track & ,
G4double ,
G4ForceCondition * condition )
virtual

Implements G4VProcess.

Definition at line 216 of file G4ParallelGeometriesLimiterProcess.cc.

219{
220 // -- push previous step limitation flags and volumes:
221 // -- consider switching pointers insteads of making copies of std::vector's:
222 fParallelWorldWasLimiting = fParallelWorldIsLimiting;
223 fPreviousVolumes = fCurrentVolumes;
224
225 // -- update volumes:
226 std::size_t i = 0;
227 for ( auto navigatorIndex : fParallelWorldNavigatorIndeces )
228 {
229 fCurrentVolumes[i++] = fPathFinder->GetLocatedVolume( navigatorIndex );
230 }
231
233 return DBL_MAX;
234}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced

◆ RemoveParallelWorld()

void G4ParallelGeometriesLimiterProcess::RemoveParallelWorld ( const G4String & parallelWorldName)

Definition at line 116 of file G4ParallelGeometriesLimiterProcess.cc.

118{
119 // -- Refuse refuse removing parallel geometry during tracking time:
120 if (fIsTrackingTime)
121 {
123 ed << "G4ParallelGeometriesLimiterProcess `" << GetProcessName()
124 << "': removing a parallel world volume at tracking time is not allowed."
125 << G4endl;
126 G4Exception("G4ParallelGeometriesLimiterProcess::RemoveParallelWorld(..)",
127 "BIAS.GEN.25", JustWarning, ed, "Call ignored.");
128 return;
129 }
130 else
131 {
132 G4VPhysicalVolume* newWorld = fTransportationManager->IsWorldExisting( parallelWorldName );
133 if (newWorld == nullptr)
134 {
136 ed << "G4ParallelGeometriesLimiterProcess `" << GetProcessName()
137 << "': trying to remove an inexisting parallel world '"
138 << parallelWorldName << "'." << G4endl;
139 G4Exception("G4ParallelGeometriesLimiterProcess::RemoveParallelWorld(..)",
140 "BIAS.GEN.26", JustWarning, ed, "Call ignored.");
141 return;
142 }
143
144 // -- get position of world volume in list:
145 std::size_t iWorld = 0;
146 for ( auto knownWorld : fParallelWorlds )
147 {
148 if ( knownWorld == newWorld ) break;
149 ++iWorld;
150 }
151
152 if ( iWorld == fParallelWorlds.size() )
153 {
155 ed << "G4ParallelGeometriesLimiterProcess `" << GetProcessName()
156 << "': trying to remove an non-registerered parallel world '"
157 << parallelWorldName << "'." << G4endl;
158 G4Exception("G4ParallelGeometriesLimiterProcess::RemoveParallelWorld(..)",
159 "BIAS.GEN.27", JustWarning, ed, "Call ignored.");
160 return;
161 }
162 // -- remove from vector:
163 fParallelWorlds.erase( fParallelWorlds.begin() + iWorld );
164 }
165}

◆ SetProcessManager()

void G4ParallelGeometriesLimiterProcess::SetProcessManager ( const G4ProcessManager * mgr)
virtual

Reimplemented from G4VProcess.

Definition at line 361 of file G4ParallelGeometriesLimiterProcess.cc.

363{
364 G4BiasingProcessSharedData *sharedData(nullptr);
365
366 // -- initialize sharedData pointer:
367 if ( G4BiasingProcessSharedData::fSharedDataMap.Find(mgr) == G4BiasingProcessSharedData::fSharedDataMap.End() )
368 {
369 sharedData = new G4BiasingProcessSharedData( mgr );
370 G4BiasingProcessSharedData::fSharedDataMap[mgr] = sharedData;
371 }
372 else
373 {
374 sharedData = G4BiasingProcessSharedData::fSharedDataMap[mgr] ;
375 }
376
377 // -- add itself to the shared data:
378 if ( sharedData->fParallelGeometriesLimiterProcess == nullptr )
379 {
380 sharedData->fParallelGeometriesLimiterProcess = this;
381 }
382 else
383 {
385 ed << " Trying to add more than one G4ParallelGeometriesLimiterProcess process to the process manager "
386 << mgr
387 << " (process manager for `" << mgr->GetParticleType()->GetParticleName()
388 << "'). Only one is needed. Call ignored." << G4endl;
389 G4Exception(" G4ParallelGeometriesLimiterProcess::SetProcessManager(..)",
390 "BIAS.GEN.29", JustWarning, ed);
391 }
392}
const G4String & GetParticleName() const
G4ParticleDefinition * GetParticleType() const

◆ StartTracking()

void G4ParallelGeometriesLimiterProcess::StartTracking ( G4Track * track)
virtual

Reimplemented from G4VProcess.

Definition at line 170 of file G4ParallelGeometriesLimiterProcess.cc.

171{
172 fIsTrackingTime = true;
173
174 // -- fetch the navigators, their indeces, and activate:
175 fParallelWorldNavigators.clear();
176 fParallelWorldNavigatorIndeces.clear();
177 fParallelWorldSafeties.clear();
178 fParallelWorldIsLimiting.clear();
179 fParallelWorldWasLimiting.clear();
180 fCurrentVolumes.clear();
181 fPreviousVolumes.clear();
182 for ( auto parallelWorld : fParallelWorlds )
183 {
184 fParallelWorldNavigators.push_back( fTransportationManager->GetNavigator( parallelWorld ) );
185 fParallelWorldNavigatorIndeces.push_back( fTransportationManager->ActivateNavigator( fParallelWorldNavigators.back() ) );
186 fParallelWorldSafeties.push_back( 0.0 );
187 fParallelWorldIsLimiting.push_back( false );
188 fParallelWorldWasLimiting.push_back( false );
189 }
190
191 fPathFinder->PrepareNewTrack( track->GetPosition(), track->GetMomentumDirection() );
192 // -- Does it work at this level, after "PrepareNewTrack" above ?
193 for ( auto navigatorIndex : fParallelWorldNavigatorIndeces )
194 {
195 fPreviousVolumes.push_back( nullptr );
196 fCurrentVolumes .push_back( fPathFinder->GetLocatedVolume( navigatorIndex ) );
197 }
198
199 // -- will force updating safety:
200 fParallelWorldSafety = 0.0;
201 for ( std::size_t i = 0 ; i < fParallelWorldNavigatorIndeces.size() ; ++i )
202 {
203 fParallelWorldSafeties[i] = 0.0;
204 }
205}
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const

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