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

#include <G4ParallelWorldScoringProcess.hh>

+ Inheritance diagram for G4ParallelWorldScoringProcess:

Public Member Functions

 G4ParallelWorldScoringProcess (const G4String &processName="ParaWorldScore", G4ProcessType theType=fParameterisation)
 
virtual ~G4ParallelWorldScoringProcess ()
 
void SetParallelWorld (G4String parallelWorldName)
 
void SetParallelWorld (G4VPhysicalVolume *parallelWorld)
 
G4bool IsAtRestRequired (G4ParticleDefinition *partDef)
 
void StartTracking (G4Track *)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
void Verbose (const G4Step &) const
 
- 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
 

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
 

Detailed Description

Definition at line 68 of file G4ParallelWorldScoringProcess.hh.

Constructor & Destructor Documentation

◆ G4ParallelWorldScoringProcess()

G4ParallelWorldScoringProcess::G4ParallelWorldScoringProcess ( const G4String processName = "ParaWorldScore",
G4ProcessType  theType = fParameterisation 
)

Definition at line 53 of file G4ParallelWorldScoringProcess.cc.

55:G4VProcess(processName,theType), fGhostNavigator(0), fNavigatorID(-1), fFieldTrack('0')
56{
57 pParticleChange = &aDummyParticleChange;
58
59 fGhostStep = new G4Step();
60 fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
61 fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
62
64 fPathFinder = G4PathFinder::GetInstance();
65
66 if (verboseLevel>0)
67 {
68 G4cout << GetProcessName() << " is created " << G4endl;
69 }
70}
#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

◆ ~G4ParallelWorldScoringProcess()

G4ParallelWorldScoringProcess::~G4ParallelWorldScoringProcess ( )
virtual

Definition at line 75 of file G4ParallelWorldScoringProcess.cc.

76{
77 delete fGhostStep;
78}

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Definition at line 390 of file G4ParallelWorldScoringProcess.cc.

392{
393 // Dummy ParticleChange ie: does nothing
394 // Expecting G4Transportation to move the track
396 return pParticleChange;
397}
virtual void Initialize(const G4Track &)

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 324 of file G4ParallelWorldScoringProcess.cc.

327{
328 static G4FieldTrack endTrack('0');
329 static ELimited eLimited;
330
331 *selection = NotCandidateForSelection;
332 G4double returnedStep = DBL_MAX;
333
334 if (previousStepSize > 0.)
335 { fGhostSafety -= previousStepSize; }
336// else
337// { fGhostSafety = -1.; }
338 if (fGhostSafety < 0.) fGhostSafety = 0.0;
339
340 // ------------------------------------------
341 // Determination of the proposed STEP LENGTH:
342 // ------------------------------------------
343 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
344 {
345 // I have no chance to limit
346 returnedStep = currentMinimumStep;
347 fOnBoundary = false;
348 proposedSafety = fGhostSafety - currentMinimumStep;
349 }
350 else // (currentMinimumStep > fGhostSafety: I may limit the Step)
351 {
352 G4FieldTrackUpdator::Update(&fFieldTrack,&track);
353//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
354// ComputeStep
355//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
356 returnedStep
357 = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
358 track.GetCurrentStepNumber(),fGhostSafety,eLimited,
359 endTrack,track.GetVolume());
360//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
361 if(eLimited == kDoNot)
362 {
363 // Track is not on the boundary
364 fOnBoundary = false;
365 fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());
366 }
367 else
368 {
369 // Track is on the boundary
370 fOnBoundary = true;
371 // proposedSafety = fGhostSafety;
372 }
373 proposedSafety = fGhostSafety;
374 if(eLimited == kUnique || eLimited == kSharedOther) {
375 *selection = CandidateForSelection;
376 }else if (eLimited == kSharedTransport) {
377 returnedStep *= (1.0 + 1.0e-9);
378 // Expand to disable its selection in Step Manager comparison
379 }
380 }
381
382 // ----------------------------------------------
383 // Returns the fGhostSafety as the proposedSafety
384 // The SteppingManager will take care of keeping
385 // the smallest one.
386 // ----------------------------------------------
387 return returnedStep;
388}
@ 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 * G4ParallelWorldScoringProcess::AtRestDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4VProcess.

Definition at line 198 of file G4ParallelWorldScoringProcess.cc.

201{
202 fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
203 G4VSensitiveDetector* aSD = 0;
204 if(fOldGhostTouchable->GetVolume())
205 { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
206 fOnBoundary = false;
207 CopyStep(step);
208 fGhostPreStepPoint->SetSensitiveDetector(aSD);
209
210 fNewGhostTouchable = fOldGhostTouchable;
211
212 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
213 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
214 if(fNewGhostTouchable->GetVolume())
215 {
216 fGhostPostStepPoint->SetSensitiveDetector(
217 fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
218 }
219 else
220 { fGhostPostStepPoint->SetSensitiveDetector(0); }
221
222 if (verboseLevel>1) Verbose(step);
223
224 G4VSensitiveDetector* sd = fGhostPreStepPoint->GetSensitiveDetector();
225 if(sd)
226 {
227 sd->Hit(fGhostStep);
228 }
229
231 return pParticleChange;
232}
G4VSensitiveDetector * GetSensitiveDetector() const
void SetSensitiveDetector(G4VSensitiveDetector *)
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4TouchableHandle & GetTouchableHandle() const
G4VSensitiveDetector * GetSensitiveDetector() const
G4LogicalVolume * GetLogicalVolume() const
G4bool Hit(G4Step *aStep)
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44

◆ AtRestGetPhysicalInteractionLength()

G4double G4ParallelWorldScoringProcess::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 185 of file G4ParallelWorldScoringProcess.cc.

188{
189 *condition = Forced;
190 return DBL_MAX;
191}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced

◆ IsAtRestRequired()

G4bool G4ParallelWorldScoringProcess::IsAtRestRequired ( G4ParticleDefinition partDef)

Definition at line 109 of file G4ParallelWorldScoringProcess.cc.

111{
112 G4int pdgCode = partDef->GetPDGEncoding();
113 if(pdgCode==0)
114 {
115 G4String partName = partDef->GetParticleName();
116 if(partName=="opticalphoton") return false;
117 if(partName=="geantino") return false;
118 if(partName=="chargedgeantino") return false;
119 }
120 else
121 {
122 if(pdgCode==22) return false; // gamma
123 if(pdgCode==11) return false; // electron
124 if(pdgCode==2212) return false; // proton
125 if(pdgCode==-12) return false; // anti_nu_e
126 if(pdgCode==12) return false; // nu_e
127 if(pdgCode==-14) return false; // anti_nu_mu
128 if(pdgCode==14) return false; // nu_mu
129 if(pdgCode==-16) return false; // anti_nu_tau
130 if(pdgCode==16) return false; // nu_tau
131 }
132 return true;
133}
int G4int
Definition: G4Types.hh:66
const G4String & GetParticleName() const

Referenced by G4RunManager::ConstructScoringWorlds().

◆ PostStepDoIt()

G4VParticleChange * G4ParallelWorldScoringProcess::PostStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4VProcess.

Definition at line 255 of file G4ParallelWorldScoringProcess.cc.

258{
259 fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
260 G4VSensitiveDetector* aSD = 0;
261 if(fOldGhostTouchable->GetVolume())
262 { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
263 CopyStep(step);
264 fGhostPreStepPoint->SetSensitiveDetector(aSD);
265
266 // fPathFinder->Locate( track.GetPosition(),
267 // track.GetMomentumDirection(),
268 // true);
269
270 // fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
271 // step.GetPostStepPoint()->GetMomentumDirection());
272
273 if(fOnBoundary)
274 {
275//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
276// Locate the point and get new touchable
277//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
278 //?? fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
279 //?? step.GetPostStepPoint()->GetMomentumDirection());
280 fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
281//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
282 }
283 else
284 {
285// Do I need this ??????????????????????????????????????????????????????????
286// fGhostNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
287// ?????????????????????????????????????????????????????????????????????????
288
289 // fPathFinder->ReLocate(track.GetPosition());
290
291 // reuse the touchable
292 fNewGhostTouchable = fOldGhostTouchable;
293 }
294
295 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
296 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
297
298 if(fNewGhostTouchable->GetVolume())
299 {
300 fGhostPostStepPoint->SetSensitiveDetector(
301 fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
302 }
303 else
304 { fGhostPostStepPoint->SetSensitiveDetector(0); }
305
306 if (verboseLevel>1) Verbose(step);
307
308 G4VSensitiveDetector* sd = fGhostPreStepPoint->GetSensitiveDetector();
309 if(sd)
310 {
311 sd->Hit(fGhostStep);
312 }
313
314 pParticleChange->Initialize(track); // Does not change the track properties
315 return pParticleChange;
316}
G4TouchableHandle CreateTouchableHandle(G4int navId) const

◆ PostStepGetPhysicalInteractionLength()

G4double G4ParallelWorldScoringProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 240 of file G4ParallelWorldScoringProcess.cc.

244{
245 // I must be invoked anyway to score the hit.
247 return DBL_MAX;
248}
@ StronglyForced

◆ SetParallelWorld() [1/2]

void G4ParallelWorldScoringProcess::SetParallelWorld ( G4String  parallelWorldName)

Definition at line 85 of file G4ParallelWorldScoringProcess.cc.

87{
88//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
89// Get pointers of the parallel world and its navigator
90//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
91 fGhostWorldName = parallelWorldName;
92 fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
93 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
94//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
95}
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
G4Navigator * GetNavigator(const G4String &worldName)

Referenced by G4RunManager::ConstructScoringWorlds().

◆ SetParallelWorld() [2/2]

void G4ParallelWorldScoringProcess::SetParallelWorld ( G4VPhysicalVolume parallelWorld)

Definition at line 97 of file G4ParallelWorldScoringProcess.cc.

99{
100//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
101// Get pointer of navigator
102//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
103 fGhostWorldName = parallelWorld->GetName();
104 fGhostWorld = parallelWorld;
105 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
106//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
107}
const G4String & GetName() const

◆ StartTracking()

void G4ParallelWorldScoringProcess::StartTracking ( G4Track trk)
virtual

Reimplemented from G4VProcess.

Definition at line 141 of file G4ParallelWorldScoringProcess.cc.

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

◆ Verbose()

void G4ParallelWorldScoringProcess::Verbose ( const G4Step step) const

Definition at line 424 of file G4ParallelWorldScoringProcess.cc.

425{
426 G4cout << "In mass geometry ------------------------------------------------" << G4endl;
427 G4cout << " StepLength : " << step.GetStepLength()/mm << " TotalEnergyDeposit : "
428 << step.GetTotalEnergyDeposit()/MeV << G4endl;
429 G4cout << " PreStepPoint : "
430 << step.GetPreStepPoint()->GetPhysicalVolume()->GetName() << " - ";
433 else
434 { G4cout << "NoProcessAssigned"; }
435 G4cout << G4endl;
436 G4cout << " " << step.GetPreStepPoint()->GetPosition() << G4endl;
437 G4cout << " PostStepPoint : ";
440 else
441 { G4cout << "OutOfWorld"; }
442 G4cout << " - ";
445 else
446 { G4cout << "NoProcessAssigned"; }
447 G4cout << G4endl;
448 G4cout << " " << step.GetPostStepPoint()->GetPosition() << G4endl;
449
450 G4cout << "In ghost geometry ------------------------------------------------" << G4endl;
451 G4cout << " StepLength : " << fGhostStep->GetStepLength()/mm
452 << " TotalEnergyDeposit : "
453 << fGhostStep->GetTotalEnergyDeposit()/MeV << G4endl;
454 G4cout << " PreStepPoint : "
455 << fGhostStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() << " ["
456 << fGhostStep->GetPreStepPoint()->GetTouchable()->GetReplicaNumber()
457 << " ]" << " - ";
458 if(fGhostStep->GetPreStepPoint()->GetProcessDefinedStep())
460 else
461 { G4cout << "NoProcessAssigned"; }
462 G4cout << G4endl;
463 G4cout << " " << fGhostStep->GetPreStepPoint()->GetPosition() << G4endl;
464 G4cout << " PostStepPoint : ";
465 if(fGhostStep->GetPostStepPoint()->GetPhysicalVolume())
466 {
467 G4cout << fGhostStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() << " ["
468 << fGhostStep->GetPostStepPoint()->GetTouchable()->GetReplicaNumber()
469 << " ]";
470 }
471 else
472 { G4cout << "OutOfWorld"; }
473 G4cout << " - ";
474 if(fGhostStep->GetPostStepPoint()->GetProcessDefinedStep())
476 else
477 { G4cout << "NoProcessAssigned"; }
478 G4cout << G4endl;
479 G4cout << " " << fGhostStep->GetPostStepPoint()->GetPosition() << " == "
480 << fGhostStep->GetTrack()->GetMomentumDirection()
481 << G4endl;
482
483}
const G4VTouchable * GetTouchable() const
const G4VProcess * GetProcessDefinedStep() const
const G4ThreeVector & GetPosition() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4Track * GetTrack() const
G4double GetStepLength() const
G4double GetTotalEnergyDeposit() const
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:58

Referenced by AtRestDoIt(), and PostStepDoIt().


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