Geant4  10.02.p03
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 G4VParticleChange * PostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChange * AtRestDoIt (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 G4VParticleChange * AlongStepDoIt (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 ()
 

Private Member Functions

 XWrapperDiscreteProcess (XWrapperDiscreteProcess &)
 
XWrapperDiscreteProcessoperator= (const XWrapperDiscreteProcess &)
 

Private Attributes

G4int bBothOrCrystalOrDetectorPhysics
 Decide whether to use nuclei (+1) or electron (-1) or both (0) More...
 
G4int bNucleiOrElectronFlag
 
G4VDiscreteProcessfRegisteredProcess
 
G4ParticleChangeForNothingfParticleChangeForNothing
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
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 38 of file XWrapperDiscreteProcess.hh.

Constructor & Destructor Documentation

◆ XWrapperDiscreteProcess() [1/3]

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

Definition at line 42 of file XWrapperDiscreteProcess.cc.

43 :G4VDiscreteProcess(aName){
44  if (verboseLevel>1) {
45  G4cout << GetProcessName() << " is created "<< G4endl;
46  }
50 }
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
G4int bBothOrCrystalOrDetectorPhysics
Decide whether to use nuclei (+1) or electron (-1) or both (0)
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForNothing * fParticleChangeForNothing
Here is the call graph for this function:
Here is the caller graph for this function:

◆ XWrapperDiscreteProcess() [2/3]

XWrapperDiscreteProcess::XWrapperDiscreteProcess ( const G4String ,
G4VDiscreteProcess toRegister 
)

Definition at line 55 of file XWrapperDiscreteProcess.cc.

57 :G4VDiscreteProcess("XWrapperDiscreteProcess"){
58  fRegisteredProcess = toRegister;
67 }
G4VDiscreteProcess * fRegisteredProcess
G4ProcessType theProcessType
Definition: G4VProcess.hh:340
G4int theProcessSubType
Definition: G4VProcess.hh:343
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
G4bool isPostStepDoItIsEnabled() const
Definition: G4VProcess.hh:532
G4bool isAlongStepDoItIsEnabled() const
Definition: G4VProcess.hh:526
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
G4int bBothOrCrystalOrDetectorPhysics
Decide whether to use nuclei (+1) or electron (-1) or both (0)
G4bool isAtRestDoItIsEnabled() const
Definition: G4VProcess.hh:520
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:414
G4ParticleChangeForNothing * fParticleChangeForNothing
Here is the call graph for this function:

◆ ~XWrapperDiscreteProcess()

XWrapperDiscreteProcess::~XWrapperDiscreteProcess ( )
virtual

Definition at line 71 of file XWrapperDiscreteProcess.cc.

71  {
72  ;
73 }
Here is the call graph for this function:

◆ XWrapperDiscreteProcess() [3/3]

XWrapperDiscreteProcess::XWrapperDiscreteProcess ( XWrapperDiscreteProcess right)
private

Definition at line 78 of file XWrapperDiscreteProcess.cc.

78  :
79 G4VDiscreteProcess(right){
80  ;
81 }

Member Function Documentation

◆ AtRestDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 81 of file XWrapperDiscreteProcess.hh.

82  {
83  return fRegisteredProcess->AtRestDoIt(aTrack,aStep);
84  }
G4VDiscreteProcess * fRegisteredProcess
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
Here is the call graph for this function:

◆ AtRestGetPhysicalInteractionLength()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 89 of file XWrapperDiscreteProcess.hh.

90  {
92  condition);
93  };
G4double condition(const G4ErrorSymMatrix &m)
G4VDiscreteProcess * fRegisteredProcess
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
Here is the call graph for this function:

◆ BuildPhysicsTable()

void XWrapperDiscreteProcess::BuildPhysicsTable ( const G4ParticleDefinition aParticleDefinition)
virtual

Reimplemented from G4VProcess.

Definition at line 314 of file XWrapperDiscreteProcess.cc.

314  {
315  fRegisteredProcess->BuildPhysicsTable(aParticleDefinition);
316 }
G4VDiscreteProcess * fRegisteredProcess
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:210
Here is the call graph for this function:
Here is the caller graph for this function:

◆ BuildWorkerPhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 130 of file XWrapperDiscreteProcess.hh.

G4VDiscreteProcess * fRegisteredProcess
virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition &part)
Definition: G4VProcess.cc:202
Here is the call graph for this function:

◆ DumpInfo()

virtual void XWrapperDiscreteProcess::DumpInfo ( ) const
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 123 of file XWrapperDiscreteProcess.hh.

G4VDiscreteProcess * fRegisteredProcess
virtual void DumpInfo() const
Definition: G4VProcess.cc:178
Here is the call graph for this function:

◆ EndTracking()

virtual void XWrapperDiscreteProcess::EndTracking ( )
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 114 of file XWrapperDiscreteProcess.hh.

G4VDiscreteProcess * fRegisteredProcess
virtual void EndTracking()
Definition: G4VProcess.cc:113
Here is the call graph for this function:

◆ GetDensity()

G4double XWrapperDiscreteProcess::GetDensity ( const G4Track &  aTrack)

Definition at line 149 of file XWrapperDiscreteProcess.cc.

149  {
150  //Retrieve nuclei and electron density
151  //from ExExChParticleUserInfo object
152  G4double vDensity = 1.;
153 
154  if(ItHasToWork(aTrack)){
155  ExExChParticleUserInfo* chanInfo =
156  (ExExChParticleUserInfo*) aTrack.GetUserInformation();
157 
158  if(chanInfo){
159  if(bNucleiOrElectronFlag == +1){
160  vDensity = chanInfo->GetNucleiDensity();
161  }
162  else if(bNucleiOrElectronFlag == -1){
163  vDensity = chanInfo->GetElectronDensity();
164  }
165  else{
166  vDensity = (chanInfo->GetNucleiDensity()
167  + chanInfo->GetElectronDensity())/2.;
168  }
169  }
170  }
171 
172  return vDensity;
173 }
G4int ItHasToWork(const G4Track &)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetDensityPreviousStep()

G4double XWrapperDiscreteProcess::GetDensityPreviousStep ( const G4Track &  aTrack)

Definition at line 178 of file XWrapperDiscreteProcess.cc.

178  {
179  //Retrieve nuclei and electron density
180  //from ExExChParticleUserInfo object
181 
182  G4double vDensityPreviousStep = 1.;
183 
184  if(ItHasToWork(aTrack)){
185  ExExChParticleUserInfo* chanInfo =
186  (ExExChParticleUserInfo*) aTrack.GetUserInformation();
187 
188  if(bNucleiOrElectronFlag == +1){
189  vDensityPreviousStep = chanInfo->GetNucleiDensityPreviousStep();
190  }
191  else if(bNucleiOrElectronFlag == -1){
192  vDensityPreviousStep = chanInfo->GetElectronDensityPreviousStep();
193  }
194  else{
195  vDensityPreviousStep =
196  (chanInfo->GetNucleiDensityPreviousStep()
197  + chanInfo->GetElectronDensityPreviousStep())/2.;
198  }
199  }
200 
201  return vDensityPreviousStep;
202 }
G4int ItHasToWork(const G4Track &)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 231 of file XWrapperDiscreteProcess.cc.

233  {
234  return DBL_MAX;
235 }
#define DBL_MAX
Definition: templates.hh:83
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetNucleiOrElectronFlag()

G4int XWrapperDiscreteProcess::GetNucleiOrElectronFlag ( )

Definition at line 117 of file XWrapperDiscreteProcess.cc.

117  {
118  return bNucleiOrElectronFlag;
119 }

◆ GetProcessManager()

virtual const G4ProcessManager* XWrapperDiscreteProcess::GetProcessManager ( )
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 118 of file XWrapperDiscreteProcess.hh.

G4VDiscreteProcess * fRegisteredProcess
virtual const G4ProcessManager * GetProcessManager()
Definition: G4VProcess.hh:514
Here is the call graph for this function:

◆ IsApplicable()

G4bool XWrapperDiscreteProcess::IsApplicable ( const G4ParticleDefinition aParticleDefinition)
virtual

Reimplemented from G4VProcess.

Definition at line 307 of file XWrapperDiscreteProcess.cc.

307  {
308  return fRegisteredProcess->IsApplicable(aParticleDefinition);
309 }
G4VDiscreteProcess * fRegisteredProcess
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:205
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ItHasToWork()

G4int XWrapperDiscreteProcess::ItHasToWork ( const G4Track &  aTrack)

Definition at line 123 of file XWrapperDiscreteProcess.cc.

123  {
124  ExExChParticleUserInfo* chanInfo =
125  (ExExChParticleUserInfo*) aTrack.GetUserInformation();
126 
127  if(chanInfo){
128  if((chanInfo->GetInTheCrystal() == true) &&
131  return 1;
132  }
133  if((chanInfo->GetInTheCrystal() == false) &&
136  return 2;
137  }
138  }
139  else {
140  G4cout << G4endl << "XWrapperDiscreteProcess::";
141  G4cout << "ERROR - no ExExChParticleUserInfo object Detected";
142  G4cout << G4endl;
143  }
144 
145  return 0;
146 }
G4GLOB_DLL std::ostream G4cout
G4int bBothOrCrystalOrDetectorPhysics
Decide whether to use nuclei (+1) or electron (-1) or both (0)
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

XWrapperDiscreteProcess& XWrapperDiscreteProcess::operator= ( const XWrapperDiscreteProcess )
private

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 291 of file XWrapperDiscreteProcess.cc.

292  {
293  if(ItHasToWork(aTrack) == 1){
294  return fRegisteredProcess->PostStepDoIt(aTrack, aStep);
295  }
296  else if(ItHasToWork(aTrack) == 2){
297  return fRegisteredProcess->PostStepDoIt(aTrack, aStep);
298  }
300  return pParticleChange;
301 }
G4VDiscreteProcess * fRegisteredProcess
G4int ItHasToWork(const G4Track &)
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4ParticleChangeForNothing * fParticleChangeForNothing
Here is the call graph for this function:

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 240 of file XWrapperDiscreteProcess.cc.

242  {
243  if(ItHasToWork(aTrack) == 1){
244  G4double vDensity = GetDensity(aTrack);
245  G4double vDensityPreviousStep = GetDensityPreviousStep(aTrack);
246 
247  if ( (previousStepSize < 0.0) ||
249  // beginning of tracking (or just after DoIt of this process)
251  } else if ( previousStepSize > 0.0) {
252  // subtract NumberOfInteractionLengthLeft
254  * vDensityPreviousStep);
255  } else {
256  // zero step DO NOTHING
257  }
258 
259  G4double regIntLength =
261  previousStepSize * vDensityPreviousStep,
262  condition);
263  G4double regIntNumber =
265  if(regIntNumber!=0){
266  currentInteractionLength = regIntLength / regIntNumber;
267  }
268  else{
269  return DBL_MAX;
270  }
271  theNumberOfInteractionLengthLeft = regIntNumber;
272 
275  if ( vDensity == 0. ) return DBL_MAX;
276  currentInteractionLength /= vDensity;
278  }
279  else if(ItHasToWork(aTrack) == 2){
281  previousStepSize,
282  condition);
283  }
284  else{
285  return DBL_MAX;
286  }
287 }
G4double condition(const G4ErrorSymMatrix &m)
G4VDiscreteProcess * fRegisteredProcess
G4double GetNumberOfInteractionLengthLeft() const
Definition: G4VProcess.hh:453
G4int ItHasToWork(const G4Track &)
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
G4double GetDensityPreviousStep(const G4Track &)
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:
Here is the caller graph for this function:

◆ PreparePhysicsTable()

void XWrapperDiscreteProcess::PreparePhysicsTable ( const G4ParticleDefinition aParticleDefinition)
virtual

Reimplemented from G4VProcess.

Definition at line 321 of file XWrapperDiscreteProcess.cc.

321  {
322  fRegisteredProcess->PreparePhysicsTable(aParticleDefinition);
323 }
G4VDiscreteProcess * fRegisteredProcess
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:217
Here is the call graph for this function:
Here is the caller graph for this function:

◆ PrepareWorkerPhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 132 of file XWrapperDiscreteProcess.hh.

G4VDiscreteProcess * fRegisteredProcess
virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.cc:207
Here is the call graph for this function:

◆ RegisterProcess() [1/2]

void XWrapperDiscreteProcess::RegisterProcess ( G4VDiscreteProcess toRegister)

Definition at line 85 of file XWrapperDiscreteProcess.cc.

85  {
86  fRegisteredProcess = toRegister;
92 }
G4VDiscreteProcess * fRegisteredProcess
G4ProcessType theProcessType
Definition: G4VProcess.hh:340
G4int theProcessSubType
Definition: G4VProcess.hh:343
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
G4bool isPostStepDoItIsEnabled() const
Definition: G4VProcess.hh:532
G4bool isAlongStepDoItIsEnabled() const
Definition: G4VProcess.hh:526
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
G4bool isAtRestDoItIsEnabled() const
Definition: G4VProcess.hh:520
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:414
Here is the call graph for this function:
Here is the caller graph for this function:

◆ RegisterProcess() [2/2]

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

Definition at line 96 of file XWrapperDiscreteProcess.cc.

98  {
99  fRegisteredProcess = toRegister;
100  bNucleiOrElectronFlag = flag;
107 }
G4VDiscreteProcess * fRegisteredProcess
G4ProcessType theProcessType
Definition: G4VProcess.hh:340
G4int theProcessSubType
Definition: G4VProcess.hh:343
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
G4bool isPostStepDoItIsEnabled() const
Definition: G4VProcess.hh:532
G4bool isAlongStepDoItIsEnabled() const
Definition: G4VProcess.hh:526
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
G4int bBothOrCrystalOrDetectorPhysics
Decide whether to use nuclei (+1) or electron (-1) or both (0)
G4bool isAtRestDoItIsEnabled() const
Definition: G4VProcess.hh:520
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:414
Here is the call graph for this function:

◆ ResetNumberOfInteractionLengthLeft()

virtual void XWrapperDiscreteProcess::ResetNumberOfInteractionLengthLeft ( )
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 120 of file XWrapperDiscreteProcess.hh.

G4VDiscreteProcess * fRegisteredProcess
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:95
Here is the call graph for this function:
Here is the caller graph for this function:

◆ RetrievePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 339 of file XWrapperDiscreteProcess.cc.

341  {
342  return fRegisteredProcess->RetrievePhysicsTable(aParticleDefinition,
343  aString,
344  aBool);
345 }
G4VDiscreteProcess * fRegisteredProcess
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:236
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetMasterProcess()

virtual void XWrapperDiscreteProcess::SetMasterProcess ( G4VProcess masterP)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 125 of file XWrapperDiscreteProcess.hh.

125  {
127  static_cast<XWrapperDiscreteProcess*>(masterP)->fRegisteredProcess
128  );
129  };
virtual void SetMasterProcess(G4VProcess *masterP)
Definition: G4VProcess.cc:212
G4VDiscreteProcess * fRegisteredProcess
Here is the call graph for this function:

◆ SetNucleiOrElectronFlag()

void XWrapperDiscreteProcess::SetNucleiOrElectronFlag ( G4int  flag)

Definition at line 111 of file XWrapperDiscreteProcess.cc.

111  {
112  bNucleiOrElectronFlag = flag;
113 }

◆ SetProcessManager()

virtual void XWrapperDiscreteProcess::SetProcessManager ( const G4ProcessManager aPM)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 116 of file XWrapperDiscreteProcess.hh.

G4VDiscreteProcess * fRegisteredProcess
virtual void SetProcessManager(const G4ProcessManager *)
Definition: G4VProcess.hh:508
Here is the call graph for this function:

◆ StartTracking()

void XWrapperDiscreteProcess::StartTracking ( G4Track *  aTrack)
virtual

Reimplemented from G4VProcess.

Definition at line 222 of file XWrapperDiscreteProcess.cc.

222  {
227 }
G4VDiscreteProcess * fRegisteredProcess
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:
Here is the caller graph for this function:

◆ StorePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 328 of file XWrapperDiscreteProcess.cc.

330  {
331  return fRegisteredProcess->StorePhysicsTable(aParticleDefinition,
332  aString,
333  aBool);
334 }
G4VDiscreteProcess * fRegisteredProcess
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:231
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ bBothOrCrystalOrDetectorPhysics

G4int XWrapperDiscreteProcess::bBothOrCrystalOrDetectorPhysics
private

Decide whether to use nuclei (+1) or electron (-1) or both (0)

Definition at line 68 of file XWrapperDiscreteProcess.hh.

◆ bNucleiOrElectronFlag

G4int XWrapperDiscreteProcess::bNucleiOrElectronFlag
private

Definition at line 69 of file XWrapperDiscreteProcess.hh.

◆ fParticleChangeForNothing

G4ParticleChangeForNothing* XWrapperDiscreteProcess::fParticleChangeForNothing
private

Definition at line 71 of file XWrapperDiscreteProcess.hh.

◆ fRegisteredProcess

G4VDiscreteProcess* XWrapperDiscreteProcess::fRegisteredProcess
private

Definition at line 70 of file XWrapperDiscreteProcess.hh.


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