Geant4  10.02.p03
XWrapperContinuousDiscreteProcess Class Reference

#include <XWrapperContinuousDiscreteProcess.hh>

Inheritance diagram for XWrapperContinuousDiscreteProcess:
Collaboration diagram for XWrapperContinuousDiscreteProcess:

Public Member Functions

 XWrapperContinuousDiscreteProcess (const G4String &processName="XWrapperContinuousDiscreteProcess")
 
 XWrapperContinuousDiscreteProcess (const G4String &processName, G4VContinuousDiscreteProcess *)
 
 XWrapperContinuousDiscreteProcess (const G4String &processName, G4ProcessType)
 
G4int ItHasToWork (const G4Track &)
 
virtual ~XWrapperContinuousDiscreteProcess ()
 
void RegisterProcess (G4VContinuousDiscreteProcess *)
 
void RegisterProcess (G4VContinuousDiscreteProcess *, G4int, G4int aBool=0)
 
G4VContinuousDiscreteProcessGetRegisteredProcess ()
 
G4double GetDensity (const G4Track &)
 
G4double GetDensityPreviousStep (const G4Track &)
 
void SetNucleiOrElectronFlag (G4int)
 
G4int GetNucleiOrElectronFlag ()
 
virtual G4VParticleChange * PostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChange * AlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChange * AtRestDoIt (const G4Track &aTrack, const G4Step &aStep)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &aTrack, G4ForceCondition *condition)
 
virtual void StartTracking (G4Track *aTrack)
 
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)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &aPD)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &aPD)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &aPD)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *aPD, const G4String &aString, G4bool aBool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *aPD, const G4String &aString, G4bool aBool)
 
- Public Member Functions inherited from G4VContinuousDiscreteProcess
 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
 
virtual ~G4VContinuousDiscreteProcess ()
 
- 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 *)
 
virtual G4double GetContinuousStepLimit (const G4Track &, G4double, G4double, G4double &)
 
- Protected Member Functions inherited from G4VContinuousDiscreteProcess
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Private Member Functions

 XWrapperContinuousDiscreteProcess (XWrapperContinuousDiscreteProcess &)
 
XWrapperContinuousDiscreteProcessoperator= (const XWrapperContinuousDiscreteProcess &right)
 

Private Attributes

G4int bBothOrCrystalOrDetectorPhysics
 
G4int bNucleiOrElectronFlag
 
G4VContinuousDiscreteProcessfRegisteredProcess
 
const G4Step theStepCopy
 
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 XWrapperContinuousDiscreteProcess.hh.

Constructor & Destructor Documentation

◆ XWrapperContinuousDiscreteProcess() [1/4]

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

Definition at line 44 of file XWrapperContinousDiscreteProcess.cc.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ XWrapperContinuousDiscreteProcess() [2/4]

XWrapperContinuousDiscreteProcess::XWrapperContinuousDiscreteProcess ( const G4String processName,
G4VContinuousDiscreteProcess toRegister 
)

Definition at line 57 of file XWrapperContinousDiscreteProcess.cc.

60  fRegisteredProcess = toRegister;
61  if (verboseLevel>1) {
62  G4cout << GetProcessName() << " is created "<< G4endl;
63  }
72 }
G4ProcessType theProcessType
Definition: G4VProcess.hh:340
G4int verboseLevel
Definition: G4VProcess.hh:368
G4int theProcessSubType
Definition: G4VProcess.hh:343
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
G4bool isPostStepDoItIsEnabled() const
Definition: G4VProcess.hh:532
G4bool isAlongStepDoItIsEnabled() const
Definition: G4VProcess.hh:526
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
G4ParticleChangeForNothing * fParticleChangeForNothing
G4bool isAtRestDoItIsEnabled() const
Definition: G4VProcess.hh:520
#define G4endl
Definition: G4ios.hh:61
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:414
Here is the call graph for this function:

◆ XWrapperContinuousDiscreteProcess() [3/4]

XWrapperContinuousDiscreteProcess::XWrapperContinuousDiscreteProcess ( const G4String processName,
G4ProcessType  aProcessType 
)

Definition at line 77 of file XWrapperContinousDiscreteProcess.cc.

79 :G4VContinuousDiscreteProcess(aName,aProcessType){
80  if (verboseLevel>1) {
81  G4cout << GetProcessName() << " is created "<< G4endl;
82  }
83  theProcessType = aProcessType;
87 }
G4ProcessType theProcessType
Definition: G4VProcess.hh:340
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForNothing * fParticleChangeForNothing
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

◆ ~XWrapperContinuousDiscreteProcess()

XWrapperContinuousDiscreteProcess::~XWrapperContinuousDiscreteProcess ( )
virtual

Definition at line 91 of file XWrapperContinousDiscreteProcess.cc.

91  {
92  ;
93 }
Here is the call graph for this function:

◆ XWrapperContinuousDiscreteProcess() [4/4]

XWrapperContinuousDiscreteProcess::XWrapperContinuousDiscreteProcess ( XWrapperContinuousDiscreteProcess right)
private

Definition at line 98 of file XWrapperContinousDiscreteProcess.cc.

Here is the call graph for this function:

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * XWrapperContinuousDiscreteProcess::AlongStepDoIt ( const G4Track &  aTrack,
const G4Step &  aStep 
)
virtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 373 of file XWrapperContinousDiscreteProcess.cc.

374  {
375 
376  if(ItHasToWork(aTrack) == 1){
377  G4double vDensity = GetDensity(aTrack);
378  G4double vStepLengthSaved = aStep.GetStepLength();
379  const_cast<G4Step&>(aStep).SetStepLength(
380  aStep.GetStepLength() * vDensity);
382  const_cast<G4Step&>(aStep).SetStepLength(vStepLengthSaved);
383  return pParticleChange;
384  }
385  else if(ItHasToWork(aTrack) == 2){
386  return fRegisteredProcess->AlongStepDoIt(aTrack, aStep);
387  }
389  return pParticleChange;
390 }
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
G4ParticleChangeForNothing * fParticleChangeForNothing
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ AlongStepGetPhysicalInteractionLength()

G4double XWrapperContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength ( const G4Track &  aTrack,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety,
G4GPILSelection *  selection 
)
virtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 321 of file XWrapperContinousDiscreteProcess.cc.

325  {
326 
327  if(ItHasToWork(aTrack) == 1){
328  G4double vDensityPreviousStep = GetDensityPreviousStep(aTrack);
329 
331  previousStepSize * vDensityPreviousStep,
332  currentMinimumStep,
333  currentSafety,
334  selection);
335  }
336  else if(ItHasToWork(aTrack) == 2){
338  previousStepSize,
339  currentMinimumStep,
340  currentSafety,
341  selection);
342  }
343  else{
344  return DBL_MAX;
345  }
346 }
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
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:

◆ AtRestDoIt()

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 88 of file XWrapperContinuousDiscreteProcess.hh.

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

◆ AtRestGetPhysicalInteractionLength()

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 104 of file XWrapperContinuousDiscreteProcess.hh.

105  {
107  condition);
108  };
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
Here is the call graph for this function:

◆ BuildPhysicsTable()

virtual void XWrapperContinuousDiscreteProcess::BuildPhysicsTable ( const G4ParticleDefinition aPD)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 137 of file XWrapperContinuousDiscreteProcess.hh.

virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:210
Here is the call graph for this function:

◆ BuildWorkerPhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 128 of file XWrapperContinuousDiscreteProcess.hh.

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

◆ DumpInfo()

virtual void XWrapperContinuousDiscreteProcess::DumpInfo ( ) const
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 119 of file XWrapperContinuousDiscreteProcess.hh.

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

◆ EndTracking()

virtual void XWrapperContinuousDiscreteProcess::EndTracking ( )
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 112 of file XWrapperContinuousDiscreteProcess.hh.

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

◆ GetContinuousStepLimit()

G4double XWrapperContinuousDiscreteProcess::GetContinuousStepLimit ( const G4Track &  ,
G4double  ,
G4double  ,
G4double  
)
protectedvirtual

Implements G4VContinuousDiscreteProcess.

Definition at line 395 of file XWrapperContinousDiscreteProcess.cc.

398  {
399  return DBL_MAX;
400 }
#define DBL_MAX
Definition: templates.hh:83
Here is the caller graph for this function:

◆ GetDensity()

G4double XWrapperContinuousDiscreteProcess::GetDensity ( const G4Track &  aTrack)

Definition at line 171 of file XWrapperContinousDiscreteProcess.cc.

171  {
172  //Retrieve nuclei and electron density
173  //from ExExChParticleUserInfo object
174  G4double vDensity = 1.;
175 
176  if(ItHasToWork(aTrack) == 1){
177  ExExChParticleUserInfo* chanInfo =
178  (ExExChParticleUserInfo*) aTrack.GetUserInformation();
179 
180  if(chanInfo){
181  if(bNucleiOrElectronFlag == +1){
182  vDensity = chanInfo->GetNucleiDensity();
183  }
184  else if(bNucleiOrElectronFlag == -1){
185  vDensity = chanInfo->GetElectronDensity();
186  }
187  else{
188  vDensity = (chanInfo->GetNucleiDensity()
189  + chanInfo->GetElectronDensity())/2.;
190  }
191  }
192  }
193 
194  return vDensity;
195 }
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetDensityPreviousStep()

G4double XWrapperContinuousDiscreteProcess::GetDensityPreviousStep ( const G4Track &  aTrack)

Definition at line 200 of file XWrapperContinousDiscreteProcess.cc.

200  {
201  //Retrieve nuclei and electron density
202  //from ExExChParticleUserInfo object
203 
204  G4double vDensityPreviousStep = 1.;
205 
206  if(ItHasToWork(aTrack) == 1){
207  ExExChParticleUserInfo* chanInfo =
208  (ExExChParticleUserInfo*) aTrack.GetUserInformation();
209 
210  if(bNucleiOrElectronFlag == +1){
211  vDensityPreviousStep = chanInfo->GetNucleiDensityPreviousStep();
212  }
213  else if(bNucleiOrElectronFlag == -1){
214  vDensityPreviousStep = chanInfo->GetElectronDensityPreviousStep();
215  }
216  else{
217  vDensityPreviousStep =
218  (chanInfo->GetNucleiDensityPreviousStep()
219  + chanInfo->GetElectronDensityPreviousStep())/2.;
220  }
221  }
222 
223  return vDensityPreviousStep;
224 }
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetMeanFreePath()

G4double XWrapperContinuousDiscreteProcess::GetMeanFreePath ( const G4Track &  aTrack,
G4double  previousStepSize,
G4ForceCondition *  condition 
)
protectedvirtual

Implements G4VContinuousDiscreteProcess.

Definition at line 254 of file XWrapperContinousDiscreteProcess.cc.

256  {
257  if(ItHasToWork(aTrack) == 2){
259  previousStepSize,
260  condition);
261  }
262  return DBL_MAX;
263 }
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
#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 XWrapperContinuousDiscreteProcess::GetNucleiOrElectronFlag ( )

Definition at line 139 of file XWrapperContinousDiscreteProcess.cc.

Here is the caller graph for this function:

◆ GetProcessManager()

virtual const G4ProcessManager* XWrapperContinuousDiscreteProcess::GetProcessManager ( )
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 115 of file XWrapperContinuousDiscreteProcess.hh.

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

◆ GetRegisteredProcess()

G4VContinuousDiscreteProcess* XWrapperContinuousDiscreteProcess::GetRegisteredProcess ( )
inline

Definition at line 57 of file XWrapperContinuousDiscreteProcess.hh.

58  {return fRegisteredProcess;};
Here is the call graph for this function:

◆ IsApplicable()

virtual G4bool XWrapperContinuousDiscreteProcess::IsApplicable ( const G4ParticleDefinition aPD)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 135 of file XWrapperContinuousDiscreteProcess.hh.

136  {return fRegisteredProcess->IsApplicable(aPD);};
virtual G4bool IsApplicable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:205
Here is the call graph for this function:

◆ ItHasToWork()

G4int XWrapperContinuousDiscreteProcess::ItHasToWork ( const G4Track &  aTrack)

Definition at line 145 of file XWrapperContinousDiscreteProcess.cc.

145  {
146  ExExChParticleUserInfo* chanInfo =
147  (ExExChParticleUserInfo*) aTrack.GetUserInformation();
148 
149  if(chanInfo){
150  if((chanInfo->GetInTheCrystal() == true) &&
153  return 1;
154  }
155  if((chanInfo->GetInTheCrystal() == false) &&
158  return 2;
159  }
160  }
161  else {
162  G4cout << G4endl << "XWrapperDiscreteProcess::";
163  G4cout << "ERROR - no ExExChParticleUserInfo object Detected";
164  G4cout << G4endl;
165  }
166 
167  return 0;
168 }
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:

◆ operator=()

XWrapperContinuousDiscreteProcess& XWrapperContinuousDiscreteProcess::operator= ( const XWrapperContinuousDiscreteProcess right)
private
Here is the caller graph for this function:

◆ PostStepDoIt()

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 351 of file XWrapperContinousDiscreteProcess.cc.

352  {
353 
354  if(ItHasToWork(aTrack) == 1){
355  G4double vDensity = GetDensity(aTrack);
356  G4double vStepLengthSaved = aStep.GetStepLength();
357  const_cast<G4Step&>(aStep).SetStepLength(
358  aStep.GetStepLength() * vDensity);
360  const_cast<G4Step&>(aStep).SetStepLength(vStepLengthSaved);
361  return pParticleChange;
362  }
363  else if(ItHasToWork(aTrack) == 2){
364  return fRegisteredProcess->PostStepDoIt(aTrack, aStep);
365  }
367  return pParticleChange;
368 }
G4ParticleChangeForNothing * fParticleChangeForNothing
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 268 of file XWrapperContinousDiscreteProcess.cc.

270  {
271 
272  if(ItHasToWork(aTrack) == 1){
273  G4double vDensity = GetDensity(aTrack);
274  G4double vDensityPreviousStep = GetDensityPreviousStep(aTrack);
275 
276  if ( (previousStepSize < 0.0) ||
278  // beginning of tracking (or just after DoIt of this process)
280  } else if ( previousStepSize > 0.0) {
281  // subtract NumberOfInteractionLengthLeft
283  * vDensityPreviousStep);
284  } else {
285  // zero step DO NOTHING
286  }
287 
288  G4double regIntLength =
290  previousStepSize * vDensityPreviousStep,
291  condition);
292  G4double regIntNumber =
294  if(regIntNumber!=0){
295  currentInteractionLength = regIntLength / regIntNumber;
296  }
297  else{
298  return DBL_MAX;
299  }
300  theNumberOfInteractionLengthLeft = regIntNumber;
301 
304  if ( vDensity == 0. ) return DBL_MAX;
305  currentInteractionLength /= vDensity;
307  }
308  else if(ItHasToWork(aTrack) == 2){
310  previousStepSize,
311  condition);
312  }
313  else{
314  return DBL_MAX;
315  }
316 }
G4double condition(const G4ErrorSymMatrix &m)
G4double GetNumberOfInteractionLengthLeft() const
Definition: G4VProcess.hh:453
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
G4double currentInteractionLength
Definition: G4VProcess.hh:297
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
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()

virtual void XWrapperContinuousDiscreteProcess::PreparePhysicsTable ( const G4ParticleDefinition aPD)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 139 of file XWrapperContinuousDiscreteProcess.hh.

virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:217
Here is the call graph for this function:

◆ PrepareWorkerPhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 130 of file XWrapperContinuousDiscreteProcess.hh.

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

◆ RegisterProcess() [1/2]

void XWrapperContinuousDiscreteProcess::RegisterProcess ( G4VContinuousDiscreteProcess toRegister)

Definition at line 106 of file XWrapperContinousDiscreteProcess.cc.

106  {
107  fRegisteredProcess = toRegister;
113 }
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 XWrapperContinuousDiscreteProcess::RegisterProcess ( G4VContinuousDiscreteProcess toRegister,
G4int  flag,
G4int  aBool = 0 
)

Definition at line 118 of file XWrapperContinousDiscreteProcess.cc.

120  {
121  fRegisteredProcess = toRegister;
127  bNucleiOrElectronFlag = flag;
129 }
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:

◆ ResetNumberOfInteractionLengthLeft()

virtual void XWrapperContinuousDiscreteProcess::ResetNumberOfInteractionLengthLeft ( )
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 117 of file XWrapperContinuousDiscreteProcess.hh.

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

◆ RetrievePhysicsTable()

virtual G4bool XWrapperContinuousDiscreteProcess::RetrievePhysicsTable ( const G4ParticleDefinition aPD,
const G4String aString,
G4bool  aBool 
)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 145 of file XWrapperContinuousDiscreteProcess.hh.

148  {return fRegisteredProcess->RetrievePhysicsTable(aPD,aString,aBool);};
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:236
Here is the call graph for this function:

◆ SetMasterProcess()

virtual void XWrapperContinuousDiscreteProcess::SetMasterProcess ( G4VProcess masterP)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 122 of file XWrapperContinuousDiscreteProcess.hh.

122  {
124  static_cast<XWrapperContinuousDiscreteProcess*>(masterP)
126  };
virtual void SetMasterProcess(G4VProcess *masterP)
Definition: G4VProcess.cc:212
Here is the call graph for this function:

◆ SetNucleiOrElectronFlag()

void XWrapperContinuousDiscreteProcess::SetNucleiOrElectronFlag ( G4int  flag)

Definition at line 133 of file XWrapperContinousDiscreteProcess.cc.

Here is the caller graph for this function:

◆ SetProcessManager()

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

Reimplemented from G4VProcess.

Definition at line 113 of file XWrapperContinuousDiscreteProcess.hh.

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

◆ StartTracking()

void XWrapperContinuousDiscreteProcess::StartTracking ( G4Track *  aTrack)
virtual

Reimplemented from G4VProcess.

Definition at line 244 of file XWrapperContinousDiscreteProcess.cc.

244  {
249 }
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()

virtual G4bool XWrapperContinuousDiscreteProcess::StorePhysicsTable ( const G4ParticleDefinition aPD,
const G4String aString,
G4bool  aBool 
)
inlinevirtual

Reimplemented from G4VProcess.

Definition at line 141 of file XWrapperContinuousDiscreteProcess.hh.

144  {return fRegisteredProcess->StorePhysicsTable(aPD,aString,aBool);};
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &, G4bool)
Definition: G4VProcess.hh:231
Here is the call graph for this function:

Member Data Documentation

◆ bBothOrCrystalOrDetectorPhysics

G4int XWrapperContinuousDiscreteProcess::bBothOrCrystalOrDetectorPhysics
private

Definition at line 73 of file XWrapperContinuousDiscreteProcess.hh.

◆ bNucleiOrElectronFlag

G4int XWrapperContinuousDiscreteProcess::bNucleiOrElectronFlag
private

Definition at line 74 of file XWrapperContinuousDiscreteProcess.hh.

◆ fParticleChangeForNothing

G4ParticleChangeForNothing* XWrapperContinuousDiscreteProcess::fParticleChangeForNothing
private

Definition at line 80 of file XWrapperContinuousDiscreteProcess.hh.

◆ fRegisteredProcess

G4VContinuousDiscreteProcess* XWrapperContinuousDiscreteProcess::fRegisteredProcess
private

Definition at line 77 of file XWrapperContinuousDiscreteProcess.hh.

◆ theStepCopy

const G4Step XWrapperContinuousDiscreteProcess::theStepCopy
private

Definition at line 79 of file XWrapperContinuousDiscreteProcess.hh.


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