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

#include <XWrapperDiscreteProcess.hh>

Inheritance diagram for XWrapperDiscreteProcess:
Collaboration diagram for XWrapperDiscreteProcess:

Public Member Functions

 XWrapperDiscreteProcess (const G4String &processName="XWrapperDiscreteProcess")
 
 XWrapperDiscreteProcess (const G4String &, G4VDiscreteProcess *)
 
G4int ItHasToWork (const G4Track &)
 
virtual ~XWrapperDiscreteProcess ()
 
void RegisterProcess (G4VDiscreteProcess *)
 
void RegisterProcess (G4VDiscreteProcess *, G4int, G4int aBool=0)
 
G4double GetDensity (const G4Track &)
 
G4double GetDensityPreviousStep (const G4Track &)
 
void SetNucleiOrElectronFlag (G4int)
 
G4int GetNucleiOrElectronFlag ()
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &aTrack, const G4Step &aStep)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &aTrack, G4ForceCondition *condition)
 
void StartTracking (G4Track *)
 
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)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *aPM)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
virtual void DumpInfo () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &aPD)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &aPD)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
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)
 
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)
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
const G4VProcessGetMasterProcess () const
 

Protected Member Functions

virtual G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- 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 41 of file XWrapperDiscreteProcess.hh.

Constructor & Destructor Documentation

XWrapperDiscreteProcess::XWrapperDiscreteProcess ( const G4String processName = "XWrapperDiscreteProcess")

Definition at line 45 of file XWrapperDiscreteProcess.cc.

46 :G4VDiscreteProcess(aName){
47  if (verboseLevel>1) {
48  G4cout << GetProcessName() << " is created "<< G4endl;
49  }
50  bNucleiOrElectronFlag = +0;
51  bBothOrCrystalOrDetectorPhysics = +0;
52  fParticleChangeForNothing = new G4ParticleChangeForNothing();
53 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4GLOB_DLL std::ostream G4cout
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

XWrapperDiscreteProcess::XWrapperDiscreteProcess ( const G4String ,
G4VDiscreteProcess toRegister 
)

Definition at line 58 of file XWrapperDiscreteProcess.cc.

60 :G4VDiscreteProcess("XWrapperDiscreteProcess"){
61  fRegisteredProcess = toRegister;
62  theProcessType = fRegisteredProcess->GetProcessType();
63  theProcessSubType = fRegisteredProcess->GetProcessSubType();
64  enableAtRestDoIt = fRegisteredProcess->isAtRestDoItIsEnabled();
65  enableAlongStepDoIt = fRegisteredProcess->isAlongStepDoItIsEnabled();
66  enablePostStepDoIt = fRegisteredProcess->isPostStepDoItIsEnabled();
67  bNucleiOrElectronFlag = +0;
68  bBothOrCrystalOrDetectorPhysics = +0;
69  fParticleChangeForNothing = new G4ParticleChangeForNothing();
70 }
G4ProcessType theProcessType
Definition: G4VProcess.hh:340
G4bool isAlongStepDoItIsEnabled() const
Definition: G4VProcess.hh:526
G4bool isPostStepDoItIsEnabled() const
Definition: G4VProcess.hh:532
G4bool isAtRestDoItIsEnabled() const
Definition: G4VProcess.hh:520
G4int theProcessSubType
Definition: G4VProcess.hh:343
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:414
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426

Here is the call graph for this function:

XWrapperDiscreteProcess::~XWrapperDiscreteProcess ( )
virtual

Definition at line 74 of file XWrapperDiscreteProcess.cc.

74  {
75  ;
76 }

Member Function Documentation

virtual G4VParticleChange* XWrapperDiscreteProcess::AtRestDoIt ( const G4Track aTrack,
const G4Step aStep 
)
inlinevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 84 of file XWrapperDiscreteProcess.hh.

85  {
86  return fRegisteredProcess->AtRestDoIt(aTrack,aStep);
87  }
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)

Here is the call graph for this function:

virtual G4double XWrapperDiscreteProcess::AtRestGetPhysicalInteractionLength ( const G4Track aTrack,
G4ForceCondition condition 
)
inlinevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 92 of file XWrapperDiscreteProcess.hh.

93  {
94  return fRegisteredProcess->AtRestGetPhysicalInteractionLength(aTrack,
95  condition);
96  };
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)

Here is the call graph for this function:

void XWrapperDiscreteProcess::BuildPhysicsTable ( const G4ParticleDefinition aParticleDefinition)
virtual

Reimplemented from G4VProcess.

Definition at line 317 of file XWrapperDiscreteProcess.cc.

317  {
318  fRegisteredProcess->BuildPhysicsTable(aParticleDefinition);
319 }
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:210

Here is the call graph for this function:

virtual void XWrapperDiscreteProcess::BuildWorkerPhysicsTable ( const G4ParticleDefinition aPD)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 133 of file XWrapperDiscreteProcess.hh.

134  {fRegisteredProcess->BuildWorkerPhysicsTable(aPD);};
virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition &part)
Definition: G4VProcess.cc:202

Here is the call graph for this function:

virtual void XWrapperDiscreteProcess::DumpInfo ( ) const
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 126 of file XWrapperDiscreteProcess.hh.

126 {fRegisteredProcess->DumpInfo();};
virtual void DumpInfo() const
Definition: G4VProcess.cc:178

Here is the call graph for this function:

virtual void XWrapperDiscreteProcess::EndTracking ( )
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 117 of file XWrapperDiscreteProcess.hh.

117 {fRegisteredProcess->EndTracking();};
virtual void EndTracking()
Definition: G4VProcess.cc:113

Here is the call graph for this function:

G4double XWrapperDiscreteProcess::GetDensity ( const G4Track aTrack)

Definition at line 152 of file XWrapperDiscreteProcess.cc.

152  {
153  //Retrieve nuclei and electron density
154  //from ExExChParticleUserInfo object
155  G4double vDensity = 1.;
156 
157  if(ItHasToWork(aTrack)){
158  ExExChParticleUserInfo* chanInfo =
160 
161  if(chanInfo){
162  if(bNucleiOrElectronFlag == +1){
163  vDensity = chanInfo->GetNucleiDensity();
164  }
165  else if(bNucleiOrElectronFlag == -1){
166  vDensity = chanInfo->GetElectronDensity();
167  }
168  else{
169  vDensity = (chanInfo->GetNucleiDensity()
170  + chanInfo->GetElectronDensity())/2.;
171  }
172  }
173  }
174 
175  return vDensity;
176 }
G4int ItHasToWork(const G4Track &)
G4VUserTrackInformation * GetUserInformation() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4double XWrapperDiscreteProcess::GetDensityPreviousStep ( const G4Track aTrack)

Definition at line 181 of file XWrapperDiscreteProcess.cc.

181  {
182  //Retrieve nuclei and electron density
183  //from ExExChParticleUserInfo object
184 
185  G4double vDensityPreviousStep = 1.;
186 
187  if(ItHasToWork(aTrack)){
188  ExExChParticleUserInfo* chanInfo =
190 
191  if(bNucleiOrElectronFlag == +1){
192  vDensityPreviousStep = chanInfo->GetNucleiDensityPreviousStep();
193  }
194  else if(bNucleiOrElectronFlag == -1){
195  vDensityPreviousStep = chanInfo->GetElectronDensityPreviousStep();
196  }
197  else{
198  vDensityPreviousStep =
199  (chanInfo->GetNucleiDensityPreviousStep()
200  + chanInfo->GetElectronDensityPreviousStep())/2.;
201  }
202  }
203 
204  return vDensityPreviousStep;
205 }
G4int ItHasToWork(const G4Track &)
G4VUserTrackInformation * GetUserInformation() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4double XWrapperDiscreteProcess::GetMeanFreePath ( const G4Track ,
G4double  ,
G4ForceCondition  
)
protectedvirtual

Implements G4VDiscreteProcess.

Definition at line 234 of file XWrapperDiscreteProcess.cc.

236  {
237  return DBL_MAX;
238 }
#define DBL_MAX
Definition: templates.hh:83
G4int XWrapperDiscreteProcess::GetNucleiOrElectronFlag ( )

Definition at line 120 of file XWrapperDiscreteProcess.cc.

120  {
121  return bNucleiOrElectronFlag;
122 }
virtual const G4ProcessManager* XWrapperDiscreteProcess::GetProcessManager ( )
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 121 of file XWrapperDiscreteProcess.hh.

122  {return fRegisteredProcess->GetProcessManager();};
virtual const G4ProcessManager * GetProcessManager()
Definition: G4VProcess.hh:514

Here is the call graph for this function:

G4bool XWrapperDiscreteProcess::IsApplicable ( const G4ParticleDefinition aParticleDefinition)
virtual

Reimplemented from G4VProcess.

Definition at line 310 of file XWrapperDiscreteProcess.cc.

310  {
311  return fRegisteredProcess->IsApplicable(aParticleDefinition);
312 }
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:205

Here is the call graph for this function:

G4int XWrapperDiscreteProcess::ItHasToWork ( const G4Track aTrack)

Definition at line 126 of file XWrapperDiscreteProcess.cc.

126  {
127  ExExChParticleUserInfo* chanInfo =
129 
130  if(chanInfo){
131  if((chanInfo->GetInTheCrystal() == true) &&
132  (bBothOrCrystalOrDetectorPhysics == 1 ||
133  bBothOrCrystalOrDetectorPhysics == 0)){
134  return 1;
135  }
136  if((chanInfo->GetInTheCrystal() == false) &&
137  (bBothOrCrystalOrDetectorPhysics == 2 ||
138  bBothOrCrystalOrDetectorPhysics == 0)){
139  return 2;
140  }
141  }
142  else {
143  G4cout << G4endl << "XWrapperDiscreteProcess::";
144  G4cout << "ERROR - no ExExChParticleUserInfo object Detected";
145  G4cout << G4endl;
146  }
147 
148  return 0;
149 }
G4VUserTrackInformation * GetUserInformation() const
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VDiscreteProcess.

Definition at line 294 of file XWrapperDiscreteProcess.cc.

295  {
296  if(ItHasToWork(aTrack) == 1){
297  return fRegisteredProcess->PostStepDoIt(aTrack, aStep);
298  }
299  else if(ItHasToWork(aTrack) == 2){
300  return fRegisteredProcess->PostStepDoIt(aTrack, aStep);
301  }
302  pParticleChange = fParticleChangeForNothing;
303  return pParticleChange;
304 }
G4int ItHasToWork(const G4Track &)
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)

Here is the call graph for this function:

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

Reimplemented from G4VDiscreteProcess.

Definition at line 243 of file XWrapperDiscreteProcess.cc.

245  {
246  if(ItHasToWork(aTrack) == 1){
247  G4double vDensity = GetDensity(aTrack);
248  G4double vDensityPreviousStep = GetDensityPreviousStep(aTrack);
249 
250  if ( (previousStepSize < 0.0) ||
252  // beginning of tracking (or just after DoIt of this process)
254  } else if ( previousStepSize > 0.0) {
255  // subtract NumberOfInteractionLengthLeft
257  * vDensityPreviousStep);
258  } else {
259  // zero step DO NOTHING
260  }
261 
262  G4double regIntLength =
263  fRegisteredProcess->PostStepGetPhysicalInteractionLength(aTrack,
264  previousStepSize * vDensityPreviousStep,
265  condition);
266  G4double regIntNumber =
267  fRegisteredProcess->GetNumberOfInteractionLengthLeft();
268  if(regIntNumber!=0){
269  currentInteractionLength = regIntLength / regIntNumber;
270  }
271  else{
272  return DBL_MAX;
273  }
274  theNumberOfInteractionLengthLeft = regIntNumber;
275 
278  if ( vDensity == 0. ) return DBL_MAX;
279  currentInteractionLength /= vDensity;
281  }
282  else if(ItHasToWork(aTrack) == 2){
283  return fRegisteredProcess->PostStepGetPhysicalInteractionLength(aTrack,
284  previousStepSize,
285  condition);
286  }
287  else{
288  return DBL_MAX;
289  }
290 }
G4double condition(const G4ErrorSymMatrix &m)
G4int ItHasToWork(const G4Track &)
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
G4double GetDensityPreviousStep(const G4Track &)
G4double GetNumberOfInteractionLengthLeft() const
Definition: G4VProcess.hh:453
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4double currentInteractionLength
Definition: G4VProcess.hh:297
G4double GetDensity(const G4Track &)
virtual void ResetNumberOfInteractionLengthLeft()
void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
Definition: G4VProcess.hh:544
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:

void XWrapperDiscreteProcess::PreparePhysicsTable ( const G4ParticleDefinition aParticleDefinition)
virtual

Reimplemented from G4VProcess.

Definition at line 324 of file XWrapperDiscreteProcess.cc.

324  {
325  fRegisteredProcess->PreparePhysicsTable(aParticleDefinition);
326 }
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:217

Here is the call graph for this function:

virtual void XWrapperDiscreteProcess::PrepareWorkerPhysicsTable ( const G4ParticleDefinition aPD)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 135 of file XWrapperDiscreteProcess.hh.

136  {fRegisteredProcess->PrepareWorkerPhysicsTable(aPD);};
virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.cc:207

Here is the call graph for this function:

void XWrapperDiscreteProcess::RegisterProcess ( G4VDiscreteProcess toRegister)

Definition at line 88 of file XWrapperDiscreteProcess.cc.

88  {
89  fRegisteredProcess = toRegister;
90  theProcessType = fRegisteredProcess->GetProcessType();
91  theProcessSubType = fRegisteredProcess->GetProcessSubType();
92  enableAtRestDoIt = fRegisteredProcess->isAtRestDoItIsEnabled();
93  enableAlongStepDoIt = fRegisteredProcess->isAlongStepDoItIsEnabled();
94  enablePostStepDoIt = fRegisteredProcess->isPostStepDoItIsEnabled();
95 }
G4ProcessType theProcessType
Definition: G4VProcess.hh:340
G4bool isAlongStepDoItIsEnabled() const
Definition: G4VProcess.hh:526
G4bool isPostStepDoItIsEnabled() const
Definition: G4VProcess.hh:532
G4bool isAtRestDoItIsEnabled() const
Definition: G4VProcess.hh:520
G4int theProcessSubType
Definition: G4VProcess.hh:343
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:414
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426

Here is the call graph for this function:

Here is the caller graph for this function:

void XWrapperDiscreteProcess::RegisterProcess ( G4VDiscreteProcess toRegister,
G4int  flag,
G4int  aBool = 0 
)

Definition at line 99 of file XWrapperDiscreteProcess.cc.

101  {
102  fRegisteredProcess = toRegister;
103  bNucleiOrElectronFlag = flag;
104  bBothOrCrystalOrDetectorPhysics = region;
105  theProcessType = fRegisteredProcess->GetProcessType();
106  theProcessSubType = fRegisteredProcess->GetProcessSubType();
107  enableAtRestDoIt = fRegisteredProcess->isAtRestDoItIsEnabled();
108  enableAlongStepDoIt = fRegisteredProcess->isAlongStepDoItIsEnabled();
109  enablePostStepDoIt = fRegisteredProcess->isPostStepDoItIsEnabled();
110 }
G4ProcessType theProcessType
Definition: G4VProcess.hh:340
G4bool isAlongStepDoItIsEnabled() const
Definition: G4VProcess.hh:526
G4bool isPostStepDoItIsEnabled() const
Definition: G4VProcess.hh:532
G4bool isAtRestDoItIsEnabled() const
Definition: G4VProcess.hh:520
G4int theProcessSubType
Definition: G4VProcess.hh:343
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:414
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426

Here is the call graph for this function:

virtual void XWrapperDiscreteProcess::ResetNumberOfInteractionLengthLeft ( )
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 123 of file XWrapperDiscreteProcess.hh.

124  {fRegisteredProcess->ResetNumberOfInteractionLengthLeft();};
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:95

Here is the call graph for this function:

Here is the caller graph for this function:

G4bool XWrapperDiscreteProcess::RetrievePhysicsTable ( const G4ParticleDefinition aParticleDefinition,
const G4String aString,
G4bool  aBool 
)
virtual

Reimplemented from G4VProcess.

Definition at line 342 of file XWrapperDiscreteProcess.cc.

344  {
345  return fRegisteredProcess->RetrievePhysicsTable(aParticleDefinition,
346  aString,
347  aBool);
348 }
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:236

Here is the call graph for this function:

virtual void XWrapperDiscreteProcess::SetMasterProcess ( G4VProcess masterP)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 128 of file XWrapperDiscreteProcess.hh.

128  {
129  fRegisteredProcess->SetMasterProcess(
130  static_cast<XWrapperDiscreteProcess*>(masterP)->fRegisteredProcess
131  );
132  };
virtual void SetMasterProcess(G4VProcess *masterP)
Definition: G4VProcess.cc:212

Here is the call graph for this function:

void XWrapperDiscreteProcess::SetNucleiOrElectronFlag ( G4int  flag)

Definition at line 114 of file XWrapperDiscreteProcess.cc.

114  {
115  bNucleiOrElectronFlag = flag;
116 }
virtual void XWrapperDiscreteProcess::SetProcessManager ( const G4ProcessManager aPM)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 119 of file XWrapperDiscreteProcess.hh.

120  {fRegisteredProcess->SetProcessManager(aPM);};
virtual void SetProcessManager(const G4ProcessManager *)
Definition: G4VProcess.hh:508

Here is the call graph for this function:

void XWrapperDiscreteProcess::StartTracking ( G4Track aTrack)
virtual

Reimplemented from G4VProcess.

Definition at line 225 of file XWrapperDiscreteProcess.cc.

225  {
226  fRegisteredProcess->StartTracking(aTrack);
230 }
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:101
G4double currentInteractionLength
Definition: G4VProcess.hh:297
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:300

Here is the call graph for this function:

G4bool XWrapperDiscreteProcess::StorePhysicsTable ( const G4ParticleDefinition aParticleDefinition,
const G4String aString,
G4bool  aBool 
)
virtual

Reimplemented from G4VProcess.

Definition at line 331 of file XWrapperDiscreteProcess.cc.

333  {
334  return fRegisteredProcess->StorePhysicsTable(aParticleDefinition,
335  aString,
336  aBool);
337 }
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:231

Here is the call graph for this function:


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