Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4WeightCutOffProcess Class Reference

#include <G4WeightCutOffProcess.hh>

Inheritance diagram for G4WeightCutOffProcess:
G4VProcess

Public Member Functions

 G4WeightCutOffProcess (G4double wsurvival, G4double wlimit, G4double isource, G4VIStore *istore, const G4String &aName="WeightCutOffProcess", G4bool para=false)
 
virtual ~G4WeightCutOffProcess ()
 
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 &)
 
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
 
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 &)
 

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 59 of file G4WeightCutOffProcess.hh.

Constructor & Destructor Documentation

G4WeightCutOffProcess::G4WeightCutOffProcess ( G4double  wsurvival,
G4double  wlimit,
G4double  isource,
G4VIStore istore,
const G4String aName = "WeightCutOffProcess",
G4bool  para = false 
)

Definition at line 55 of file G4WeightCutOffProcess.cc.

References FatalException, G4cout, G4endl, G4Exception(), G4PathFinder::GetInstance(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4VProcess::GetProcessName(), G4TransportationManager::GetTransportationManager(), G4VProcess::pParticleChange, and G4VProcess::verboseLevel.

61  : G4VProcess(aName),
62  fParticleChange(new G4ParticleChange),
63  fWeightSurvival(wsurvival),
64  fWeightLimit(wlimit),
65  fSourceImportance(isource),
66  fIStore(istore),
67  // fGCellFinder(aGCellFinder),
68  fGhostWorldName("NoParallelWorld"), fGhostWorld(0),
69  fGhostNavigator(0), fNavigatorID(-1), fFieldTrack('0')
70 {
71  if (!fParticleChange)
72  {
73  G4Exception("G4WeightCutOffProcess::G4WeightCutOffProcess()",
74  "FatalError", FatalException,
75  "Failed to allocate G4ParticleChange !");
76  }
77 
78  G4VProcess::pParticleChange = fParticleChange;
79 
80  fGhostStep = new G4Step();
81  fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
82  fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
83 
84  fTransportationManager = G4TransportationManager::GetTransportationManager();
85  fPathFinder = G4PathFinder::GetInstance();
86 
87  if (verboseLevel>0)
88  {
89  G4cout << GetProcessName() << " is created " << G4endl;
90  }
91 
92  paraflag = para;
93 
94 }
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
G4StepPoint * GetPreStepPoint() const
G4GLOB_DLL std::ostream G4cout
Definition: G4Step.hh:76
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4TransportationManager * GetTransportationManager()
G4StepPoint * GetPostStepPoint() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
#define G4endl
Definition: G4ios.hh:61
G4WeightCutOffProcess::~G4WeightCutOffProcess ( )
virtual

Definition at line 96 of file G4WeightCutOffProcess.cc.

97 {
98  delete fParticleChange;
99  // delete fGhostStep;
100 }

Member Function Documentation

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

Implements G4VProcess.

Definition at line 388 of file G4WeightCutOffProcess.cc.

References G4VParticleChange::Initialize(), and G4VProcess::pParticleChange.

389 {
390  // Dummy ParticleChange ie: does nothing
391  // Expecting G4Transportation to move the track
392  pParticleChange->Initialize(track);
393  return pParticleChange;
394 
395  // return 0;
396 }
virtual void Initialize(const G4Track &)
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4double G4WeightCutOffProcess::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
)
virtual

Implements G4VProcess.

Definition at line 299 of file G4WeightCutOffProcess.cc.

References CandidateForSelection, G4Navigator::ComputeSafety(), G4PathFinder::ComputeStep(), DBL_MAX, G4ThreadLocal, G4Track::GetCurrentStepNumber(), G4Track::GetVolume(), kDoNot, kSharedOther, kSharedTransport, kUnique, NotCandidateForSelection, and G4FieldTrackUpdator::Update().

302 {
303  if(paraflag) {
304  static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ; if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ; G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
305  static G4ThreadLocal ELimited *eLimited_G4MT_TLS_ = 0 ; if (!eLimited_G4MT_TLS_) eLimited_G4MT_TLS_ = new ELimited ; ELimited &eLimited = *eLimited_G4MT_TLS_;
306 
307  *selection = NotCandidateForSelection;
308  G4double returnedStep = DBL_MAX;
309 
310  if (previousStepSize > 0.)
311  { fGhostSafety -= previousStepSize; }
312  // else
313  // { fGhostSafety = -1.; }
314  if (fGhostSafety < 0.) fGhostSafety = 0.0;
315 
316  // ------------------------------------------
317  // Determination of the proposed STEP LENGTH:
318  // ------------------------------------------
319  if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
320  {
321  // I have no chance to limit
322  returnedStep = currentMinimumStep;
323  fOnBoundary = false;
324  proposedSafety = fGhostSafety - currentMinimumStep;
325  }
326  else // (currentMinimumStep > fGhostSafety: I may limit the Step)
327  {
328  G4FieldTrackUpdator::Update(&fFieldTrack,&track);
329  //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
330  // ComputeStep
331  //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
332  returnedStep
333  = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
334  track.GetCurrentStepNumber(),fGhostSafety,eLimited,
335  endTrack,track.GetVolume());
336  //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
337  if(eLimited == kDoNot)
338  {
339  // Track is not on the boundary
340  fOnBoundary = false;
341  fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());
342  }
343  else
344  {
345  // Track is on the boundary
346  fOnBoundary = true;
347  proposedSafety = fGhostSafety;
348  }
349  //xbug? proposedSafety = fGhostSafety;
350  if(eLimited == kUnique || eLimited == kSharedOther) {
351  *selection = CandidateForSelection;
352  }else if (eLimited == kSharedTransport) {
353  returnedStep *= (1.0 + 1.0e-9);
354  // Expand to disable its selection in Step Manager comparison
355  }
356 
357  }
358 
359  // ----------------------------------------------
360  // Returns the fGhostSafety as the proposedSafety
361  // The SteppingManager will take care of keeping
362  // the smallest one.
363  // ----------------------------------------------
364  return returnedStep;
365 
366  } else {
367  return DBL_MAX;
368  //not sensible! return -1.0;
369  }
370 
371 }
ELimited
#define G4ThreadLocal
Definition: tls.hh:52
G4int GetCurrentStepNumber() const
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
G4VPhysicalVolume * GetVolume() const
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 *)
G4VParticleChange * G4WeightCutOffProcess::AtRestDoIt ( const G4Track ,
const G4Step  
)
virtual

Implements G4VProcess.

Definition at line 382 of file G4WeightCutOffProcess.cc.

383 {
384  return 0;
385 }
G4double G4WeightCutOffProcess::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
)
virtual

Implements G4VProcess.

Definition at line 375 of file G4WeightCutOffProcess.cc.

377 {
378  return -1.0;
379 }
const G4String & G4WeightCutOffProcess::GetName ( void  ) const

Definition at line 293 of file G4WeightCutOffProcess.cc.

References G4VProcess::theProcessName.

294 {
295  return theProcessName;
296 }
G4String theProcessName
Definition: G4VProcess.hh:335
G4VParticleChange * G4WeightCutOffProcess::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

Implements G4VProcess.

Definition at line 190 of file G4WeightCutOffProcess.cc.

References G4PathFinder::CreateTouchableHandle(), fStopAndKill, G4UniformRand, G4VIStore::GetImportance(), G4StepPoint::GetPhysicalVolume(), G4Step::GetPostStepPoint(), G4VTouchable::GetReplicaNumber(), G4StepPoint::GetTouchable(), G4StepPoint::GetTouchableHandle(), G4Track::GetWeight(), G4ParticleChange::Initialize(), G4VParticleChange::ProposeTrackStatus(), G4VParticleChange::ProposeWeight(), and G4StepPoint::SetTouchableHandle().

192 {
193  fParticleChange->Initialize(aTrack);
194 
195  if(paraflag) {
196  fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
197  //xbug? fOnBoundary = false;
198  CopyStep(aStep);
199 
200  if(fOnBoundary)
201  {
202 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
203 // Locate the point and get new touchable
204 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
205  //?? fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
206  //?? step.GetPostStepPoint()->GetMomentumDirection());
207  fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
208 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
209  }
210  else
211  {
212  // Do I need this ??????????????????????????????????????????????????????????
213  // fGhostNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
214  // ?????????????????????????????????????????????????????????????????????????
215 
216  // fPathFinder->ReLocate(track.GetPosition());
217 
218  // reuse the touchable
219  fNewGhostTouchable = fOldGhostTouchable;
220  }
221 
222  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
223  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
224 
225  }
226 
227  if(paraflag) {
228  G4GeometryCell postCell(*(fGhostPostStepPoint->GetPhysicalVolume()),
229  fGhostPostStepPoint->GetTouchable()->GetReplicaNumber());
230 
231 
232  // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(aStep);
233  // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(fGhostStep);
234  G4double R = fSourceImportance;
235  if (fIStore)
236  {
237  G4double i = fIStore->GetImportance(postCell);
238  if (i>0)
239  {
240  R/=i;
241  }
242  }
243  G4double w = aTrack.GetWeight();
244  if (w<R*fWeightLimit)
245  {
246  G4double ws = fWeightSurvival*R;
247  G4double p = w/(ws);
248  if (G4UniformRand()<p)
249  {
250  fParticleChange->ProposeTrackStatus(fStopAndKill);
251  }
252  else
253  {
254  fParticleChange->ProposeWeight(ws);
255  }
256  }
257  } else {
258 
259  G4GeometryCell postCell(*(aStep.GetPostStepPoint()->GetPhysicalVolume()),
261 
262  // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(aStep);
263  // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(fGhostStep);
264  G4double R = fSourceImportance;
265  if (fIStore)
266  {
267  G4double i = fIStore->GetImportance(postCell);
268  if (i>0)
269  {
270  R/=i;
271  }
272  }
273  G4double w = aTrack.GetWeight();
274  if (w<R*fWeightLimit)
275  {
276  G4double ws = fWeightSurvival*R;
277  G4double p = w/(ws);
278  if (G4UniformRand()<p)
279  {
280  fParticleChange->ProposeTrackStatus(fStopAndKill);
281  }
282  else
283  {
284  fParticleChange->ProposeWeight(ws);
285  }
286  }
287  }
288 
289  return fParticleChange;
290 
291 }
const char * p
Definition: xmltok.h:285
G4TouchableHandle CreateTouchableHandle(G4int navId) const
const G4VTouchable * GetTouchable() const
void ProposeWeight(G4double finalWeight)
#define G4UniformRand()
Definition: Randomize.hh:87
G4VPhysicalVolume * GetPhysicalVolume() const
virtual G4double GetImportance(const G4GeometryCell &gCell) const =0
virtual void Initialize(const G4Track &)
G4StepPoint * GetPostStepPoint() const
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:58
G4double GetWeight() const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4TouchableHandle & GetTouchableHandle() const
G4double G4WeightCutOffProcess::PostStepGetPhysicalInteractionLength ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 178 of file G4WeightCutOffProcess.cc.

References DBL_MAX, and Forced.

180 {
181 // *condition = Forced;
182 // return kInfinity;
183 
184 // *condition = StronglyForced;
185  *condition = Forced;
186  return DBL_MAX;
187 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
void G4WeightCutOffProcess::SetParallelWorld ( G4String  parallelWorldName)

Definition at line 109 of file G4WeightCutOffProcess.cc.

References G4TransportationManager::GetNavigator(), and G4TransportationManager::GetParallelWorld().

Referenced by G4WeightCutOffConfigurator::Configure().

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

Definition at line 121 of file G4WeightCutOffProcess.cc.

References G4VPhysicalVolume::GetName(), and G4TransportationManager::GetNavigator().

122 {
123 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
124 // Get pointer of navigator
125 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
126  fGhostWorldName = parallelWorld->GetName();
127  fGhostWorld = parallelWorld;
128  fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
129 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
130 }
const G4String & GetName() const
G4Navigator * GetNavigator(const G4String &worldName)
void G4WeightCutOffProcess::StartTracking ( G4Track trk)
virtual

Reimplemented from G4VProcess.

Definition at line 137 of file G4WeightCutOffProcess.cc.

References G4TransportationManager::ActivateNavigator(), G4PathFinder::CreateTouchableHandle(), FatalException, G4Exception(), G4Track::GetMomentumDirection(), G4Track::GetPosition(), G4PathFinder::PrepareNewTrack(), and G4StepPoint::SetTouchableHandle().

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

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