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

#include <G4ImportanceProcess.hh>

+ Inheritance diagram for G4ImportanceProcess:

Public Member Functions

 G4ImportanceProcess (const G4VImportanceAlgorithm &aImportanceAlgorithm, const G4VIStore &aIstore, const G4VTrackTerminator *TrackTerminator, const G4String &aName="ImportanceProcess", G4bool para=false)
 
virtual ~G4ImportanceProcess ()
 
void SetParallelWorld (const G4String &parallelWorldName)
 
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
 
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 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 ()
 

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 57 of file G4ImportanceProcess.hh.

Constructor & Destructor Documentation

◆ G4ImportanceProcess()

G4ImportanceProcess::G4ImportanceProcess ( const G4VImportanceAlgorithm & aImportanceAlgorithm,
const G4VIStore & aIstore,
const G4VTrackTerminator * TrackTerminator,
const G4String & aName = "ImportanceProcess",
G4bool para = false )

Definition at line 49 of file G4ImportanceProcess.cc.

54 : G4VProcess(aName, fParallel),
55 fParticleChange(new G4ParticleChange),
56 fImportanceAlgorithm(aImportanceAlgorithm),
57 fIStore(aIstore),
58 fParaflag(para)
59{
60 G4cout << "### G4ImportanceProcess:: Creating " << G4endl;
61 if (TrackTerminator != nullptr)
62 {
63 fPostStepAction = new G4SamplingPostStepAction(*TrackTerminator);
64 }
65 else
66 {
67 fPostStepAction = new G4SamplingPostStepAction(*this);
68 }
69 if (fParticleChange == nullptr)
70 {
71 G4Exception("G4ImportanceProcess::G4ImportanceProcess()",
72 "FatalError", FatalException,
73 "Failed allocation of G4ParticleChange !");
74 }
75 G4VProcess::pParticleChange = fParticleChange;
76
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 G4cout << "G4ImportanceProcess:: importance process paraflag is: " << fParaflag << G4endl;
91
92}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
@ fParallel
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4PathFinder * GetInstance()
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
static G4TransportationManager * GetTransportationManager()
G4int verboseLevel
G4VParticleChange * pParticleChange
const G4String & GetProcessName() const

◆ ~G4ImportanceProcess()

G4ImportanceProcess::~G4ImportanceProcess ( )
virtual

Definition at line 94 of file G4ImportanceProcess.cc.

95{
96
97 delete fPostStepAction;
98 delete fParticleChange;
99 // delete fGhostStep;
100 // delete fGhostWorld;
101 // delete fGhostNavigator;
102
103}

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4ImportanceProcess::AlongStepDoIt ( const G4Track & aTrack,
const G4Step &  )
virtual

Implements G4VProcess.

Definition at line 440 of file G4ImportanceProcess.cc.

442{
443 // Dummy ParticleChange ie: does nothing
444 // Expecting G4Transportation to move the track
445 //AH G4cout << " along step do it " << G4endl;
447 return pParticleChange;
448}
virtual void Initialize(const G4Track &)

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 351 of file G4ImportanceProcess.cc.

355{
356 if(fParaflag) {
357 *selection = NotCandidateForSelection;
358 G4double returnedStep = DBL_MAX;
359
360 if (previousStepSize > 0.)
361 { fGhostSafety -= previousStepSize; }
362 // else
363 // { fGhostSafety = -1.; }
364 if (fGhostSafety < 0.) fGhostSafety = 0.0;
365
366 // ------------------------------------------
367 // Determination of the proposed STEP LENGTH:
368 // ------------------------------------------
369 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
370 {
371 // I have no chance to limit
372 returnedStep = currentMinimumStep;
373 fOnBoundary = false;
374 proposedSafety = fGhostSafety - currentMinimumStep;
375 //AH G4cout << " step not limited, why? " << G4endl;
376 }
377 else // (currentMinimumStep > fGhostSafety: I may limit the Step)
378 {
379 G4FieldTrackUpdator::Update(&fFieldTrack,&track);
380 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
381 // ComputeStep
382 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
383 returnedStep
384 = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
385 track.GetCurrentStepNumber(),fGhostSafety,feLimited,
386 fEndTrack,track.GetVolume());
387 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
388 if(feLimited == kDoNot)
389 {
390 //AH G4cout << " computestep came back with not-boundary " << G4endl;
391 // Track is not on the boundary
392 fOnBoundary = false;
393 fGhostSafety = fGhostNavigator->ComputeSafety(fEndTrack.GetPosition());
394 }
395 else
396 {
397 // Track is on the boundary
398 //AH G4cout << " FOUND A BOUNDARY ! " << track.GetCurrentStepNumber() << G4endl;
399 fOnBoundary = true;
400 // proposedSafety = fGhostSafety;
401 }
402 proposedSafety = fGhostSafety;
403 if(feLimited == kUnique || feLimited == kSharedOther) {
404 *selection = CandidateForSelection;
405 }else if (feLimited == kSharedTransport) {
406 returnedStep *= (1.0 + 1.0e-9);
407 // Expand to disable its selection in Step Manager comparison
408 }
409
410 }
411
412 // ----------------------------------------------
413 // Returns the fGhostSafety as the proposedSafety
414 // The SteppingManager will take care of keeping
415 // the smallest one.
416 // ----------------------------------------------
417 return returnedStep;
418
419 } else {
420
421 return DBL_MAX;
422
423 }
424
425}
@ 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 * G4ImportanceProcess::AtRestDoIt ( const G4Track & ,
const G4Step &  )
virtual

Implements G4VProcess.

Definition at line 434 of file G4ImportanceProcess.cc.

436{
437 return nullptr;
438}

◆ AtRestGetPhysicalInteractionLength()

G4double G4ImportanceProcess::AtRestGetPhysicalInteractionLength ( const G4Track & ,
G4ForceCondition *  )
virtual

Implements G4VProcess.

Definition at line 427 of file G4ImportanceProcess.cc.

430{
431 return -1.0;
432}

◆ GetName()

const G4String & G4ImportanceProcess::GetName ( ) const
virtual

Implements G4VTrackTerminator.

Definition at line 346 of file G4ImportanceProcess.cc.

347{
348 return theProcessName;
349}
G4String theProcessName

◆ KillTrack()

void G4ImportanceProcess::KillTrack ( ) const
virtual

Implements G4VTrackTerminator.

Definition at line 341 of file G4ImportanceProcess.cc.

342{
343 fParticleChange->ProposeTrackStatus(fStopAndKill);
344}
@ fStopAndKill
void ProposeTrackStatus(G4TrackStatus status)

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 202 of file G4ImportanceProcess.cc.

204{
205 fParticleChange->Initialize(aTrack);
206
207 if(aTrack.GetNextVolume() == nullptr) {
208 return fParticleChange;
209 }
210
211 if(fParaflag) {
212 fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
213 //xbug? fOnBoundary = false;
214 CopyStep(aStep);
215
216 if(fOnBoundary)
217 {
218//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
219// Locate the point and get new touchable
220//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
221 //?? fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
222 //?? step.GetPostStepPoint()->GetMomentumDirection());
223 fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
224 //AH G4cout << " on boundary " << G4endl;
225//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
226 }
227 else
228 {
229 // Do I need this ??????????????????????????????????????????????????????????
230 // fGhostNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
231 // ?????????????????????????????????????????????????????????????????????????
232
233 // fPathFinder->ReLocate(track.GetPosition());
234
235 // reuse the touchable
236 fNewGhostTouchable = fOldGhostTouchable;
237 //AH G4cout << " NOT on boundary " << G4endl;
238 }
239
240 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
241 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
242
243// if ( (aStep.GetPostStepPoint()->GetStepStatus() == fGeomBoundary)
244// && (aStep.GetStepLength() > kCarTolerance) )
245// {
246// if (aTrack.GetTrackStatus()==fStopAndKill)
247// {
248// G4cout << "WARNING - G4ImportanceProcess::PostStepDoIt()"
249// << " StopAndKill track." << G4endl;
250// }
251
252// G4StepPoint *prepoint = aStep.GetPreStepPoint();
253// G4StepPoint *postpoint = aStep.GetPostStepPoint();
254// G4GeometryCell prekey(*(prepoint->GetPhysicalVolume()),
255// prepoint->GetTouchable()->GetReplicaNumber());
256// G4GeometryCell postkey(*(postpoint->GetPhysicalVolume()),
257// postpoint->GetTouchable()->GetReplicaNumber());
258
259
260 if ( (fGhostPostStepPoint->GetStepStatus() == fGeomBoundary)
261 && (aStep.GetStepLength() > kCarTolerance) )
262 {
263 if (aTrack.GetTrackStatus()==fStopAndKill)
264 {
265 G4cout << "WARNING - G4ImportanceProcess::PostStepDoIt()"
266 << " StopAndKill track. on boundary" << G4endl;
267 }
268
269 G4GeometryCell prekey(*(fGhostPreStepPoint->GetPhysicalVolume()),
270 fGhostPreStepPoint->GetTouchable()->GetReplicaNumber());
271 G4GeometryCell postkey(*(fGhostPostStepPoint->GetPhysicalVolume()),
272 fGhostPostStepPoint->GetTouchable()->GetReplicaNumber());
273
274 //AH
275 /*
276 G4cout << G4endl;
277 G4cout << G4endl;
278 G4cout << " inside parallel importance process " << aTrack.GetCurrentStepNumber() << G4endl;
279 G4cout << G4endl;
280 G4cout << G4endl;
281 G4cout << " prekey: " << fGhostPreStepPoint->GetPhysicalVolume()->GetName() << " replica:"
282 << fGhostPreStepPoint->GetTouchable()->GetReplicaNumber() << G4endl;
283 G4cout << " prekey ISTORE: " << fIStore.GetImportance(prekey) << G4endl;
284 G4cout << " postkey: " << G4endl;
285 G4cout << " postkey ISTORE: " << fIStore.GetImportance(postkey) << G4endl;
286 */
287 //AH
288 G4Nsplit_Weight nw = fImportanceAlgorithm.
289 Calculate(fIStore.GetImportance(prekey),
290 fIStore.GetImportance(postkey),
291 aTrack.GetWeight());
292 //AH
293 /*
294 G4cout << " prekey weight: " << fIStore.GetImportance(prekey)
295 << " postkey weight: " << fIStore.GetImportance(postkey)
296 << " split weight: " << nw << G4endl;
297 G4cout << " before poststepaction " << G4endl;
298 */
299 //AH
300 fPostStepAction->DoIt(aTrack, fParticleChange, nw);
301 //AH G4cout << " after post step do it " << G4endl;
302 }
303 } else {
304 if ( (aStep.GetPostStepPoint()->GetStepStatus() == fGeomBoundary)
305 && (aStep.GetStepLength() > kCarTolerance) )
306 {
307 //AH G4cout << " inside non-parallel importance process " << G4endl;
308 if (aTrack.GetTrackStatus()==fStopAndKill)
309 {
310 G4cout << "WARNING - G4ImportanceProcess::PostStepDoIt()"
311 << " StopAndKill track. on boundary non-parallel" << G4endl;
312 }
313
314 G4StepPoint *prepoint = aStep.GetPreStepPoint();
315 G4StepPoint *postpoint = aStep.GetPostStepPoint();
316
317 G4GeometryCell prekey(*(prepoint->GetPhysicalVolume()),
318 prepoint->GetTouchable()->GetReplicaNumber());
319 G4GeometryCell postkey(*(postpoint->GetPhysicalVolume()),
320 postpoint->GetTouchable()->GetReplicaNumber());
321
322 G4Nsplit_Weight nw = fImportanceAlgorithm.
323 Calculate(fIStore.GetImportance(prekey),
324 fIStore.GetImportance(postkey),
325 aTrack.GetWeight());
326 //AH
327 /*
328 G4cout << " prekey weight: " << fIStore.GetImportance(prekey)
329 << " postkey weight: " << fIStore.GetImportance(postkey)
330 << " split weight: " << nw << G4endl;
331 G4cout << " before poststepaction 2 " << G4endl;
332 */
333 //AH
334 fPostStepAction->DoIt(aTrack, fParticleChange, nw);
335 //AH G4cout << " after poststepaction 2 " << G4endl;
336 }
337 }
338 return fParticleChange;
339}
@ fGeomBoundary
void Initialize(const G4Track &) override
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
virtual G4int GetReplicaNumber(G4int depth=0) const
G4TrackStatus GetTrackStatus() const
G4VPhysicalVolume * GetNextVolume() const
G4double GetWeight() const
virtual G4double GetImportance(const G4GeometryCell &gCell) const =0

◆ PostStepGetPhysicalInteractionLength()

G4double G4ImportanceProcess::PostStepGetPhysicalInteractionLength ( const G4Track & aTrack,
G4double previousStepSize,
G4ForceCondition * condition )
virtual

Implements G4VProcess.

Definition at line 188 of file G4ImportanceProcess.cc.

192{
193// *condition = Forced;
194// return kInfinity;
195
196// *condition = StronglyForced;
197 *condition = Forced;
198 return DBL_MAX;
199}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced

◆ SetParallelWorld()

void G4ImportanceProcess::SetParallelWorld ( const G4String & parallelWorldName)

Definition at line 112 of file G4ImportanceProcess.cc.

114{
115 G4cout << G4endl << G4endl << G4endl;
116 G4cout << "G4ImportanceProcess:: SetParallelWorld name = " << parallelWorldName << G4endl;
117//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
118// Get pointers of the parallel world and its navigator
119//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
120 fGhostWorldName = parallelWorldName;
121 fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
122 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
123//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
124}
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
G4Navigator * GetNavigator(const G4String &worldName)

Referenced by G4ImportanceConfigurator::Configure().

◆ StartTracking()

void G4ImportanceProcess::StartTracking ( G4Track * trk)
virtual

Reimplemented from G4VProcess.

Definition at line 147 of file G4ImportanceProcess.cc.

148{
149//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
150// Activate navigator and get the navigator ID
151//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
152// G4cout << " G4ParallelWorldScoringProcess::StartTracking" << G4endl;
153
154 if(fParaflag) {
155 if(fGhostNavigator != nullptr)
156 { fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator); }
157 else
158 {
159 G4Exception("G4ImportanceProcess::StartTracking",
160 "ProcParaWorld000",FatalException,
161 "G4ImportanceProcess is used for tracking without having a parallel world assigned");
162 }
163//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
164
165// G4cout << "G4ParallelWorldScoringProcess::StartTracking <<<<<<<<<<<<<<<<<< " << G4endl;
166//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
167// Let PathFinder initialize
168//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
169 fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
170//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
171
172//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
173// Setup initial touchables for the first step
174//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
175 fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
176 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
177 fNewGhostTouchable = fOldGhostTouchable;
178 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
179
180 // Initialize
181 fGhostSafety = -1.;
182 fOnBoundary = false;
183 }
184
185}
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: