Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ParallelWorldProcess Class Reference

#include <G4ParallelWorldProcess.hh>

Inheritance diagram for G4ParallelWorldProcess:
Collaboration diagram for G4ParallelWorldProcess:

Public Member Functions

 G4ParallelWorldProcess (const G4String &processName="ParaWorld", G4ProcessType theType=fParallel)
 
virtual ~G4ParallelWorldProcess ()
 
void SetParallelWorld (G4String parallelWorldName)
 
void SetParallelWorld (G4VPhysicalVolume *parallelWorld)
 
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 &, G4double, G4ForceCondition *)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
void SetLayeredMaterialFlag (G4bool flg=true)
 
G4bool GetLayeredMaterialFlag () const
 
G4bool IsAtRestRequired (G4ParticleDefinition *)
 
- 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
 
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 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
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Static Public Member Functions

static const G4StepGetHyperStep ()
 
static G4int GetHypNavigatorID ()
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Additional Inherited Members

- 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 73 of file G4ParallelWorldProcess.hh.

Constructor & Destructor Documentation

G4ParallelWorldProcess::G4ParallelWorldProcess ( const G4String processName = "ParaWorld",
G4ProcessType  theType = fParallel 
)

Definition at line 62 of file G4ParallelWorldProcess.cc.

63 :G4VProcess(processName,theType),fGhostWorld(nullptr),fGhostNavigator(nullptr),
64  fNavigatorID(-1),fFieldTrack('0'),fGhostSafety(0.),fOnBoundary(false),
65  layeredMaterialFlag(false)
66 {
67  SetProcessSubType(491);
68  if(!fpHyperStep) fpHyperStep = new G4Step();
69  iParallelWorld = ++nParallelWorlds;
70 
71  pParticleChange = &aDummyParticleChange;
72 
73  fGhostStep = new G4Step();
74  fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
75  fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
76 
77  fTransportationManager = G4TransportationManager::GetTransportationManager();
78  fTransportationManager->GetNavigatorForTracking()->SetPushVerbosity(false);
79  fPathFinder = G4PathFinder::GetInstance();
80 
81  fGhostWorldName = "** NotDefined **";
83 
84  if (verboseLevel>0)
85  {
86  G4cout << GetProcessName() << " is created " << G4endl;
87  }
88 }
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:57
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:52
G4Navigator * GetNavigatorForTracking() const
static G4ParallelWorldProcessStore * GetInstance()
G4StepPoint * GetPreStepPoint() const
G4GLOB_DLL std::ostream G4cout
void SetParallelWorld(G4ParallelWorldProcess *proc, G4String parallelWorldName)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
Definition: G4Step.hh:76
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
static G4TransportationManager * GetTransportationManager()
G4StepPoint * GetPostStepPoint() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
void SetPushVerbosity(G4bool mode)
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4ParallelWorldProcess::~G4ParallelWorldProcess ( )
virtual

Definition at line 90 of file G4ParallelWorldProcess.cc.

91 {
92  delete fGhostStep;
93  nParallelWorlds--;
94  if(nParallelWorlds==0)
95  {
96  delete fpHyperStep;
97  fpHyperStep = 0;
98  }
99 }

Member Function Documentation

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

Implements G4VProcess.

Definition at line 338 of file G4ParallelWorldProcess.cc.

340 {
341  pParticleChange->Initialize(track);
342  return pParticleChange;
343 }
virtual void Initialize(const G4Track &)
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

Here is the call graph for this function:

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

Implements G4VProcess.

Definition at line 265 of file G4ParallelWorldProcess.cc.

268 {
269  static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ; if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ; G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
270  //static ELimited eLimited;
271  ELimited eLimited;
272  ELimited eLim = kUndefLimited;
273 
274  *selection = NotCandidateForSelection;
275  G4double returnedStep = DBL_MAX;
276 
277  if (previousStepSize > 0.)
278  { fGhostSafety -= previousStepSize; }
279  if (fGhostSafety < 0.) fGhostSafety = 0.0;
280 
281  if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
282  {
283  // I have no chance to limit
284  returnedStep = currentMinimumStep;
285  fOnBoundary = false;
286  proposedSafety = fGhostSafety - currentMinimumStep;
287  eLim = kDoNot;
288  }
289  else
290  {
291  G4FieldTrackUpdator::Update(&fFieldTrack,&track);
292 
293 #ifdef G4DEBUG_PARALLEL_WORLD_PROCESS
294  if( verboseLevel > 0 ){
295  int localVerb = verboseLevel-1;
296 
297  if( localVerb == 1 ) {
298  G4cout << " Pll Wrl proc::AlongStepGPIL " << this->GetProcessName() << G4endl;
299  }else if( localVerb > 1 ) {
300  G4cout << "----------------------------------------------" << G4endl;
301  G4cout << " ParallelWorldProcess: field Track set to : " << G4endl;
302  G4cout << "----------------------------------------------" << G4endl;
303  G4cout << fFieldTrack << G4endl;
304  G4cout << "----------------------------------------------" << G4endl;
305  }
306  }
307 #endif
308 
309  returnedStep
310  = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
311  track.GetCurrentStepNumber(),fGhostSafety,eLimited,
312  endTrack,track.GetVolume());
313  if(eLimited == kDoNot)
314  {
315  fOnBoundary = false;
316  fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());
317  }
318  else
319  {
320  fOnBoundary = true;
321  // fGhostSafetyEnd = 0.0; // At end-point of expected step only
322  }
323  proposedSafety = fGhostSafety;
324  if(eLimited == kUnique || eLimited == kSharedOther) {
325  *selection = CandidateForSelection;
326  }
327  else if (eLimited == kSharedTransport) {
328  returnedStep *= (1.0 + 1.0e-9);
329  }
330  eLim = eLimited;
331  }
332 
333  if(iParallelWorld==nParallelWorlds) fNavIDHyp = 0;
334  if(eLim == kUnique || eLim == kSharedOther) fNavIDHyp = fNavigatorID;
335  return returnedStep;
336 }
G4int verboseLevel
Definition: G4VProcess.hh:368
ELimited
#define G4ThreadLocal
Definition: tls.hh:89
G4GLOB_DLL std::ostream G4cout
G4int GetCurrentStepNumber() const
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4VPhysicalVolume * GetVolume() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
#define DBL_MAX
Definition: templates.hh:83
static void Update(G4FieldTrack *, const G4Track *)

Here is the call graph for this function:

G4VParticleChange * G4ParallelWorldProcess::AtRestDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4VProcess.

Definition at line 171 of file G4ParallelWorldProcess.cc.

174 {
175 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
176 // At Rest must be registered ONLY for the particle which has other At Rest
177 // process(es).
178 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
179  fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
180  G4VSensitiveDetector* aSD = 0;
181  if(fOldGhostTouchable->GetVolume())
182  { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
183  fOnBoundary = false;
184  if(aSD)
185  {
186  CopyStep(step);
187  fGhostPreStepPoint->SetSensitiveDetector(aSD);
188 
189  fNewGhostTouchable = fOldGhostTouchable;
190 
191  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
192  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
193  if(fNewGhostTouchable->GetVolume())
194  {
195  fGhostPostStepPoint->SetSensitiveDetector(
196  fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
197  }
198  else
199  { fGhostPostStepPoint->SetSensitiveDetector(0); }
200 
201  aSD->Hit(fGhostStep);
202  }
203 
204  pParticleChange->Initialize(track);
205  return pParticleChange;
206 }
virtual void Initialize(const G4Track &)
void SetSensitiveDetector(G4VSensitiveDetector *)
G4bool Hit(G4Step *aStep)
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4VSensitiveDetector * GetSensitiveDetector() const
const G4TouchableHandle & GetTouchableHandle() const

Here is the call graph for this function:

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

Implements G4VProcess.

Definition at line 159 of file G4ParallelWorldProcess.cc.

162 {
163 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
164 // At Rest must be registered ONLY for the particle which has other At Rest
165 // process(es).
166 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
167  *condition = Forced;
168  return DBL_MAX;
169 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
const G4Step * G4ParallelWorldProcess::GetHyperStep ( )
static

Definition at line 56 of file G4ParallelWorldProcess.cc.

57 { return fpHyperStep; }

Here is the caller graph for this function:

G4int G4ParallelWorldProcess::GetHypNavigatorID ( )
static

Definition at line 58 of file G4ParallelWorldProcess.cc.

59 { return fNavIDHyp; }

Here is the caller graph for this function:

G4bool G4ParallelWorldProcess::GetLayeredMaterialFlag ( ) const
inline

Definition at line 127 of file G4ParallelWorldProcess.hh.

128  { return layeredMaterialFlag; }
G4bool G4ParallelWorldProcess::IsAtRestRequired ( G4ParticleDefinition partDef)

Definition at line 425 of file G4ParallelWorldProcess.cc.

426 {
427  G4int pdgCode = partDef->GetPDGEncoding();
428  if(pdgCode==0)
429  {
430  G4String partName = partDef->GetParticleName();
431  if(partName=="opticalphoton") return false;
432  if(partName=="geantino") return false;
433  if(partName=="chargedgeantino") return false;
434  }
435  else
436  {
437  if(pdgCode==22) return false; // gamma
438  if(pdgCode==11) return false; // electron
439  if(pdgCode==2212) return false; // proton
440  if(pdgCode==-12) return false; // anti_nu_e
441  if(pdgCode==12) return false; // nu_e
442  if(pdgCode==-14) return false; // anti_nu_mu
443  if(pdgCode==14) return false; // nu_mu
444  if(pdgCode==-16) return false; // anti_nu_tau
445  if(pdgCode==16) return false; // nu_tau
446  }
447  return true;
448 }
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Implements G4VProcess.

Definition at line 218 of file G4ParallelWorldProcess.cc.

221 {
222  fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
223  G4VSensitiveDetector* aSD = 0;
224  if(fOldGhostTouchable->GetVolume())
225  { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
226  CopyStep(step);
227  fGhostPreStepPoint->SetSensitiveDetector(aSD);
228 
229  if(fOnBoundary)
230  {
231  fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
232  }
233  else
234  {
235  fNewGhostTouchable = fOldGhostTouchable;
236  }
237 
238  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
239  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
240 
241  if(fNewGhostTouchable->GetVolume())
242  {
243  fGhostPostStepPoint->SetSensitiveDetector(
244  fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
245  }
246  else
247  { fGhostPostStepPoint->SetSensitiveDetector(0); }
248 
249  G4VSensitiveDetector* sd = fGhostPreStepPoint->GetSensitiveDetector();
250  if(sd)
251  {
252  sd->Hit(fGhostStep);
253  }
254 
255  pParticleChange->Initialize(track);
256  if(layeredMaterialFlag)
257  {
258  G4StepPoint* realWorldPostStepPoint =
259  ((G4Step*)(track.GetStep()))->GetPostStepPoint();
260  SwitchMaterial(realWorldPostStepPoint);
261  }
262  return pParticleChange;
263 }
virtual void Initialize(const G4Track &)
G4TouchableHandle CreateTouchableHandle(G4int navId) const
const G4Step * GetStep() const
void SetSensitiveDetector(G4VSensitiveDetector *)
G4bool Hit(G4Step *aStep)
Definition: G4Step.hh:76
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4VSensitiveDetector * GetSensitiveDetector() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4VSensitiveDetector * GetSensitiveDetector() const
const G4TouchableHandle & GetTouchableHandle() const

Here is the call graph for this function:

G4double G4ParallelWorldProcess::PostStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 209 of file G4ParallelWorldProcess.cc.

213 {
215  return DBL_MAX;
216 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
void G4ParallelWorldProcess::SetLayeredMaterialFlag ( G4bool  flg = true)
inline

Definition at line 125 of file G4ParallelWorldProcess.hh.

126  { layeredMaterialFlag = flg; }

Here is the caller graph for this function:

void G4ParallelWorldProcess::SetParallelWorld ( G4String  parallelWorldName)

Definition at line 102 of file G4ParallelWorldProcess.cc.

103 {
104  fGhostWorldName = parallelWorldName;
105  fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
106  fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
107  fGhostNavigator->SetPushVerbosity(false);
108 }
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
G4Navigator * GetNavigator(const G4String &worldName)
void SetPushVerbosity(G4bool mode)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ParallelWorldProcess::SetParallelWorld ( G4VPhysicalVolume parallelWorld)

Definition at line 111 of file G4ParallelWorldProcess.cc.

112 {
113  fGhostWorldName = parallelWorld->GetName();
114  fGhostWorld = parallelWorld;
115  fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
116  fGhostNavigator->SetPushVerbosity(false);
117 }
const G4String & GetName() const
G4Navigator * GetNavigator(const G4String &worldName)
void SetPushVerbosity(G4bool mode)

Here is the call graph for this function:

void G4ParallelWorldProcess::StartTracking ( G4Track trk)
virtual

Reimplemented from G4VProcess.

Definition at line 119 of file G4ParallelWorldProcess.cc.

120 {
121  if(fGhostNavigator)
122  { fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator); }
123  else
124  {
125  G4Exception("G4ParallelWorldProcess::StartTracking",
126  "ProcParaWorld000",FatalException,
127  "G4ParallelWorldProcess is used for tracking without having a parallel world assigned");
128  }
129  fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
130 
131  fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
132  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
133  fNewGhostTouchable = fOldGhostTouchable;
134  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
135 
136  fGhostSafety = -1.;
137  fOnBoundary = false;
138  fGhostPreStepPoint->SetStepStatus(fUndefined);
139  fGhostPostStepPoint->SetStepStatus(fUndefined);
140 
141 // G4VPhysicalVolume* thePhys = fNewGhostTouchable->GetVolume();
142 // if(thePhys)
143 // {
144 // G4Material* ghostMaterial = thePhys->GetLogicalVolume()->GetMaterial();
145 // if(ghostMaterial)
146 // { G4cout << " --- Material : " << ghostMaterial->GetName() << G4endl; }
147 // }
148 
149  *(fpHyperStep->GetPostStepPoint()) = *(trk->GetStep()->GetPostStepPoint());
150  if(layeredMaterialFlag)
151  {
152  G4StepPoint* realWorldPostStepPoint = trk->GetStep()->GetPostStepPoint();
153  SwitchMaterial(realWorldPostStepPoint);
154  }
155  *(fpHyperStep->GetPreStepPoint()) = *(fpHyperStep->GetPostStepPoint());
156 }
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=0)
const G4ThreeVector & GetPosition() const
G4TouchableHandle CreateTouchableHandle(G4int navId) const
const G4Step * GetStep() const
void SetStepStatus(const G4StepStatus aValue)
G4StepPoint * GetPreStepPoint() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4int ActivateNavigator(G4Navigator *aNavigator)
const G4ThreeVector & GetMomentumDirection() const
G4StepPoint * GetPostStepPoint() const
void SetTouchableHandle(const G4TouchableHandle &apValue)

Here is the call graph for this function:


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