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

#include <G4WeightWindowProcess.hh>

+ Inheritance diagram for G4WeightWindowProcess:

Public Member Functions

 G4WeightWindowProcess (const G4VWeightWindowAlgorithm &aWeightWindowAlgorithm, const G4VWeightWindowStore &aWWStore, const G4VTrackTerminator *TrackTerminator, G4PlaceOfAction placeOfAction, const G4String &aName="WeightWindowProcess", G4bool para=false)
 
virtual ~G4WeightWindowProcess ()
 
void SetParallelWorld (const G4String &parallelWorldName)
 
void SetParallelWorld (G4VPhysicalVolume *parallelWorld)
 
void StartTracking (G4Track *)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual void KillTrack () const
 
virtual const G4StringGetName () const
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- 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
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 StartTracking (G4Track *)
 
virtual void EndTracking ()
 
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 &)
 
- Public Member Functions inherited from G4VTrackTerminator
 G4VTrackTerminator ()
 
virtual ~G4VTrackTerminator ()
 
virtual void KillTrack () const =0
 
virtual const G4StringGetName () const =0
 

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
 
- Protected Attributes inherited from G4VTrackTerminator
G4double kCarTolerance
 

Detailed Description

Definition at line 58 of file G4WeightWindowProcess.hh.

Constructor & Destructor Documentation

◆ G4WeightWindowProcess()

G4WeightWindowProcess::G4WeightWindowProcess ( const G4VWeightWindowAlgorithm aWeightWindowAlgorithm,
const G4VWeightWindowStore aWWStore,
const G4VTrackTerminator TrackTerminator,
G4PlaceOfAction  placeOfAction,
const G4String aName = "WeightWindowProcess",
G4bool  para = false 
)

Definition at line 50 of file G4WeightWindowProcess.cc.

56 : G4VProcess(aName),
57 fParticleChange(new G4ParticleChange),
58 fWeightWindowAlgorithm(aWeightWindowAlgorithm),
59 fWeightWindowStore(aWWStore),
60 fPlaceOfAction(placeOfAction)
61{
62 if (TrackTerminator != nullptr)
63 {
64 fPostStepAction = new G4SamplingPostStepAction(*TrackTerminator);
65 }
66 else
67 {
68 fPostStepAction = new G4SamplingPostStepAction(*this);
69 }
70 if (fParticleChange == nullptr)
71 {
72 G4Exception("G4WeightWindowProcess::G4WeightWindowProcess()",
73 "FatalError", FatalException,
74 "Failed allocation of G4ParticleChange !");
75 }
76 G4VProcess::pParticleChange = fParticleChange;
77
78 fGhostStep = new G4Step();
79 fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
80 fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
81
83 fPathFinder = G4PathFinder::GetInstance();
84
85 if (verboseLevel>0)
86 {
87 G4cout << GetProcessName() << " is created " << G4endl;
88 }
89
90 fParaflag = para;
91
92}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:52
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
static G4TransportationManager * GetTransportationManager()
G4int verboseLevel
Definition: G4VProcess.hh:356
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

◆ ~G4WeightWindowProcess()

G4WeightWindowProcess::~G4WeightWindowProcess ( )
virtual

Definition at line 94 of file G4WeightWindowProcess.cc.

95{
96
97 delete fPostStepAction;
98 delete fParticleChange;
99 // delete fGhostStep;
100
101}

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Definition at line 379 of file G4WeightWindowProcess.cc.

381{
382 // Dummy ParticleChange ie: does nothing
383 // Expecting G4Transportation to move the track
385 return pParticleChange;
386}
virtual void Initialize(const G4Track &)

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 292 of file G4WeightWindowProcess.cc.

296{
297 if(fParaflag)
298 {
299 *selection = NotCandidateForSelection;
300 G4double returnedStep = DBL_MAX;
301
302 if (previousStepSize > 0.)
303 { fGhostSafety -= previousStepSize; }
304 // else
305 // { fGhostSafety = -1.; }
306 if (fGhostSafety < 0.) fGhostSafety = 0.0;
307
308 // ------------------------------------------
309 // Determination of the proposed STEP LENGTH:
310 // ------------------------------------------
311 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
312 {
313 // I have no chance to limit
314 returnedStep = currentMinimumStep;
315 fOnBoundary = false;
316 proposedSafety = fGhostSafety - currentMinimumStep;
317 }
318 else // (currentMinimumStep > fGhostSafety: I may limit the Step)
319 {
320 G4FieldTrackUpdator::Update(&fFieldTrack,&track);
321 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
322 // ComputeStep
323 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
324 returnedStep
325 = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
326 track.GetCurrentStepNumber(),fGhostSafety,feLimited,
327 fEndTrack,track.GetVolume());
328 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
329 if(feLimited == kDoNot)
330 {
331 // Track is not on the boundary
332 fOnBoundary = false;
333 fGhostSafety = fGhostNavigator->ComputeSafety(fEndTrack.GetPosition());
334 // fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition(), DBL_MAX, true);
335 }
336 else
337 {
338 // Track is on the boundary
339 fOnBoundary = true;
340 proposedSafety = fGhostSafety;
341 }
342 //xbug? proposedSafety = fGhostSafety;
343 if(feLimited == kUnique || feLimited == kSharedOther) {
344 *selection = CandidateForSelection;
345 }else if (feLimited == kSharedTransport) {
346 returnedStep *= (1.0 + 1.0e-9);
347 // Expand to disable its selection in Step Manager comparison
348 }
349
350 }
351
352 // ----------------------------------------------
353 // Returns the fGhostSafety as the proposedSafety
354 // The SteppingManager will take care of keeping
355 // the smallest one.
356 // ----------------------------------------------
357 return returnedStep;
358
359 } else {
360 return DBL_MAX;
361 // not sensible - time goes backwards! return -1.0;
362 }
363
364}
@ 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

◆ AtRestDoIt()

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

Implements G4VProcess.

Definition at line 373 of file G4WeightWindowProcess.cc.

375{
376 return nullptr;
377}

◆ AtRestGetPhysicalInteractionLength()

G4double G4WeightWindowProcess::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
)
virtual

Implements G4VProcess.

Definition at line 366 of file G4WeightWindowProcess.cc.

369{
370 return -1.0;
371}

◆ GetName()

const G4String & G4WeightWindowProcess::GetName ( ) const
virtual

Implements G4VTrackTerminator.

Definition at line 287 of file G4WeightWindowProcess.cc.

288{
290}

◆ KillTrack()

void G4WeightWindowProcess::KillTrack ( ) const
virtual

Implements G4VTrackTerminator.

Definition at line 282 of file G4WeightWindowProcess.cc.

283{
284 fParticleChange->ProposeTrackStatus(fStopAndKill);
285}
@ fStopAndKill
void ProposeTrackStatus(G4TrackStatus status)

◆ PostStepDoIt()

G4VParticleChange * G4WeightWindowProcess::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

Implements G4VProcess.

Definition at line 192 of file G4WeightWindowProcess.cc.

194{
195
196 fParticleChange->Initialize(aTrack);
197
198 if(fParaflag) {
199 fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
200 //xbug? fOnBoundary = false;
201 CopyStep(aStep);
202
203 if(fOnBoundary)
204 {
205//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
206// Locate the point and get new touchable
207//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
208 //?? fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
209 //?? step.GetPostStepPoint()->GetMomentumDirection());
210 fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
211//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
212 }
213 else
214 {
215 // Do I need this ??????????????????????????????????????????????????????????
216 // fGhostNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
217 // ?????????????????????????????????????????????????????????????????????????
218
219 // fPathFinder->ReLocate(track.GetPosition());
220
221 // reuse the touchable
222 fNewGhostTouchable = fOldGhostTouchable;
223 }
224
225 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
226 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
227
228 }
229
230 if (aStep.GetStepLength() > kCarTolerance)
231 {
232// if ( ( fPlaceOfAction == onBoundaryAndCollision)
233// || ( (fPlaceOfAction == onBoundary) &&
234// (aStep.GetPostStepPoint()->GetStepStatus() == fGeomBoundary) )
235// || ( (fPlaceOfAction == onCollision) &&
236// (aStep.GetPostStepPoint()->GetStepStatus() != fGeomBoundary) ) )
237 if(fParaflag) {
238 if ( ( fPlaceOfAction == onBoundaryAndCollision)
239 || ( (fPlaceOfAction == onBoundary) &&
240 (fGhostPostStepPoint->GetStepStatus() == fGeomBoundary) )
241 || ( (fPlaceOfAction == onCollision) &&
242 (fGhostPostStepPoint->GetStepStatus() != fGeomBoundary) ) )
243 {
244
245// G4StepPoint *postpoint = aStep.GetPostStepPoint();
246
247// G4GeometryCell postCell(*(postpoint->GetPhysicalVolume()),
248// postpoint->GetTouchable()->GetReplicaNumber());
249
250 G4GeometryCell postCell(*(fGhostPostStepPoint->GetPhysicalVolume()),
251 fGhostPostStepPoint->GetTouchable()->GetReplicaNumber());
252 G4Nsplit_Weight nw =
253 fWeightWindowAlgorithm.Calculate(aTrack.GetWeight(),
254 fWeightWindowStore.GetLowerWeight(postCell,
255 aTrack.GetKineticEnergy()));
256 fPostStepAction->DoIt(aTrack, fParticleChange, nw);
257 }
258 } else {
259 if ( ( fPlaceOfAction == onBoundaryAndCollision)
260 || ( (fPlaceOfAction == onBoundary) &&
262 || ( (fPlaceOfAction == onCollision) &&
264 {
265
266 G4StepPoint *postpoint = aStep.GetPostStepPoint();
267
268 G4GeometryCell postCell(*(postpoint->GetPhysicalVolume()),
269 postpoint->GetTouchable()->GetReplicaNumber());
270
271 G4Nsplit_Weight nw =
272 fWeightWindowAlgorithm.Calculate(aTrack.GetWeight(),
273 fWeightWindowStore.GetLowerWeight(postCell,
274 aTrack.GetKineticEnergy()));
275 fPostStepAction->DoIt(aTrack, fParticleChange, nw);
276 }
277 }
278 }
279 return fParticleChange;
280}
@ onCollision
@ onBoundaryAndCollision
@ onBoundary
@ fGeomBoundary
Definition: G4StepStatus.hh:43
virtual void Initialize(const G4Track &)
G4TouchableHandle CreateTouchableHandle(G4int navId) const
void DoIt(const G4Track &aTrack, G4ParticleChange *aParticleChange, const G4Nsplit_Weight &nw)
G4StepStatus GetStepStatus() const
const G4VTouchable * GetTouchable() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4TouchableHandle & GetTouchableHandle() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4double GetStepLength() const
G4double GetWeight() const
G4double GetKineticEnergy() const
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:55
virtual G4Nsplit_Weight Calculate(G4double init_w, G4double lowerWeightBound) const =0
virtual G4double GetLowerWeight(const G4GeometryCell &gCell, G4double partEnergy) const =0

◆ PostStepGetPhysicalInteractionLength()

G4double G4WeightWindowProcess::PostStepGetPhysicalInteractionLength ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 178 of file G4WeightWindowProcess.cc.

182{
183// *condition = Forced;
184// return kInfinity;
185
186// *condition = StronglyForced;
187 *condition = Forced;
188 return DBL_MAX;
189}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced

◆ SetParallelWorld() [1/2]

void G4WeightWindowProcess::SetParallelWorld ( const G4String parallelWorldName)

Definition at line 109 of file G4WeightWindowProcess.cc.

111{
112//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
113// Get pointers of the parallel world and its navigator
114//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
115 fGhostWorldName = parallelWorldName;
116 fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
117 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
118//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
119}
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
G4Navigator * GetNavigator(const G4String &worldName)

Referenced by G4WeightWindowConfigurator::Configure().

◆ SetParallelWorld() [2/2]

void G4WeightWindowProcess::SetParallelWorld ( G4VPhysicalVolume parallelWorld)

Definition at line 121 of file G4WeightWindowProcess.cc.

123{
124//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
125// Get pointer of navigator
126//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
127 fGhostWorldName = parallelWorld->GetName();
128 fGhostWorld = parallelWorld;
129 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
130//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
131}
const G4String & GetName() const

◆ StartTracking()

void G4WeightWindowProcess::StartTracking ( G4Track trk)
virtual

Reimplemented from G4VProcess.

Definition at line 138 of file G4WeightWindowProcess.cc.

139{
140//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
141// Activate navigator and get the navigator ID
142//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
143// G4cout << " G4ParallelWorldScoringProcess::StartTracking" << G4endl;
144
145 if(fParaflag) {
146 if(fGhostNavigator != nullptr)
147 { fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator); }
148 else
149 {
150 G4Exception("G4WeightWindowProcess::StartTracking",
151 "ProcParaWorld000",FatalException,
152 "G4WeightWindowProcess is used for tracking without having a parallel world assigned");
153 }
154//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
155
156// G4cout << "G4ParallelWorldScoringProcess::StartTracking <<<<<<<<<<<<<<<<<< " << G4endl;
157//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
158// Let PathFinder initialize
159//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
160 fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
161//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
162
163//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
164// Setup initial touchables for the first step
165//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
166 fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
167 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
168 fNewGhostTouchable = fOldGhostTouchable;
169 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
170
171 // Initialize
172 fGhostSafety = -1.;
173 fOnBoundary = false;
174 }
175}
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=nullptr)
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
G4int ActivateNavigator(G4Navigator *aNavigator)

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