Geant4 9.6.0
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 (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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
- 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 previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 
- Protected Attributes inherited from G4VTrackTerminator
G4double kCarTolerance
 

Detailed Description

Definition at line 61 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 54 of file G4ImportanceProcess.cc.

59 : G4VProcess(aName),
60 fParticleChange(new G4ParticleChange),
61 fImportanceAlgorithm(aImportanceAlgorithm),
62 fIStore(aIstore),
63 fPostStepAction(0),
64 fGhostNavigator(0), fNavigatorID(-1), fFieldTrack('0'),
65 paraflag(para)
66
67{
68 if (TrackTerminator)
69 {
70 fPostStepAction = new G4SamplingPostStepAction(*TrackTerminator);
71 }
72 else
73 {
74 fPostStepAction = new G4SamplingPostStepAction(*this);
75 }
76 if (!fParticleChange)
77 {
78 G4Exception("G4ImportanceProcess::G4ImportanceProcess()",
79 "FatalError", FatalException,
80 "Failed allocation of G4ParticleChange !");
81 }
82 G4VProcess::pParticleChange = fParticleChange;
83
84
85 fGhostStep = new G4Step();
86 fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
87 fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
88
90 fPathFinder = G4PathFinder::GetInstance();
91
92 if (verboseLevel>0)
93 {
94 G4cout << GetProcessName() << " is created " << G4endl;
95 }
96
97 G4cout << " importance process paraflag is: " << paraflag << G4endl;
98
99}
@ FatalException
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:57
Definition: G4Step.hh:78
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
static G4TransportationManager * GetTransportationManager()
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ ~G4ImportanceProcess()

G4ImportanceProcess::~G4ImportanceProcess ( )
virtual

Definition at line 101 of file G4ImportanceProcess.cc.

102{
103
104 delete fPostStepAction;
105 delete fParticleChange;
106 delete fGhostStep;
107
108}

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Definition at line 441 of file G4ImportanceProcess.cc.

443{
444 // Dummy ParticleChange ie: does nothing
445 // Expecting G4Transportation to move the track
446 //AH G4cout << " along step do it " << G4endl;
448 return pParticleChange;
449 // return 0;
450}
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 346 of file G4ImportanceProcess.cc.

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

◆ AtRestDoIt()

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

Implements G4VProcess.

Definition at line 435 of file G4ImportanceProcess.cc.

437{
438 return 0;
439}

◆ AtRestGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 428 of file G4ImportanceProcess.cc.

431{
432 return -1.0;
433}

◆ GetName()

const G4String & G4ImportanceProcess::GetName ( ) const
virtual

Implements G4VTrackTerminator.

Definition at line 341 of file G4ImportanceProcess.cc.

342{
343 return theProcessName;
344}
G4String theProcessName
Definition: G4VProcess.hh:335

◆ KillTrack()

void G4ImportanceProcess::KillTrack ( ) const
virtual

Implements G4VTrackTerminator.

Definition at line 336 of file G4ImportanceProcess.cc.

337{
338 fParticleChange->ProposeTrackStatus(fStopAndKill);
339}
@ fStopAndKill
void ProposeTrackStatus(G4TrackStatus status)

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 201 of file G4ImportanceProcess.cc.

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

◆ PostStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 187 of file G4ImportanceProcess.cc.

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

◆ SetParallelWorld() [1/2]

void G4ImportanceProcess::SetParallelWorld ( G4String  parallelWorldName)

Definition at line 117 of file G4ImportanceProcess.cc.

119{
120//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
121// Get pointers of the parallel world and its navigator
122//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
123 fGhostWorldName = parallelWorldName;
124 fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
125 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
126//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
127}
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
G4Navigator * GetNavigator(const G4String &worldName)

Referenced by G4ImportanceConfigurator::Configure().

◆ SetParallelWorld() [2/2]

void G4ImportanceProcess::SetParallelWorld ( G4VPhysicalVolume parallelWorld)

Definition at line 129 of file G4ImportanceProcess.cc.

131{
132//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
133// Get pointer of navigator
134//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
135 fGhostWorldName = parallelWorld->GetName();
136 fGhostWorld = parallelWorld;
137 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
138//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
139}
const G4String & GetName() const

◆ StartTracking()

void G4ImportanceProcess::StartTracking ( G4Track trk)
virtual

Reimplemented from G4VProcess.

Definition at line 146 of file G4ImportanceProcess.cc.

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

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