G4VContinuousDiscreteProcess Class Reference

#include <G4VContinuousDiscreteProcess.hh>

Inheritance diagram for G4VContinuousDiscreteProcess:

G4VProcess G4hRDEnergyLoss G4VeLowEnergyLoss G4VEnergyLoss G4VEnergyLossProcess G4VMultipleScattering G4hImpactIonisation G4alphaIonisation G4eBremsstrahlung G4eIonisation G4ePolarizedIonisation G4hBremsstrahlung G4hhIonisation G4hIonisation G4hPairProduction G4ionIonisation G4mplIonisation G4MuBremsstrahlung G4MuIonisation G4MuPairProduction G4AdjointhMultipleScattering G4eMultipleScattering G4hMultipleScattering G4MuMultipleScattering

Public Member Functions

 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
virtual ~G4VContinuousDiscreteProcess ()
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)

Protected Member Functions

virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
virtual G4double GetContinuousStepLimit (const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
void SetGPILSelection (G4GPILSelection selection)
G4GPILSelection GetGPILSelection () const

Detailed Description

Definition at line 55 of file G4VContinuousDiscreteProcess.hh.


Constructor & Destructor Documentation

G4VContinuousDiscreteProcess::G4VContinuousDiscreteProcess ( const G4String ,
G4ProcessType  aType = fNotDefined 
)

Definition at line 55 of file G4VContinuousDiscreteProcess.cc.

References G4VProcess::enableAtRestDoIt, and G4VContinuousDiscreteProcess().

Referenced by G4VContinuousDiscreteProcess().

00056   : G4VProcess(aName, aType),
00057     valueGPILSelection(CandidateForSelection)
00058 {
00059   enableAtRestDoIt = false;
00060 }

G4VContinuousDiscreteProcess::G4VContinuousDiscreteProcess ( G4VContinuousDiscreteProcess  ) 

Definition at line 66 of file G4VContinuousDiscreteProcess.cc.

References G4VContinuousDiscreteProcess().

00067                   : G4VProcess(right),
00068                     valueGPILSelection(right.valueGPILSelection)
00069 {
00070 }

G4VContinuousDiscreteProcess::~G4VContinuousDiscreteProcess (  )  [virtual]

Definition at line 62 of file G4VContinuousDiscreteProcess.cc.

00063 {
00064 }


Member Function Documentation

G4VParticleChange * G4VContinuousDiscreteProcess::AlongStepDoIt ( const G4Track ,
const G4Step  
) [virtual]

Implements G4VProcess.

Reimplemented in G4VeLowEnergyLoss, G4hImpactIonisation, G4VEnergyLoss, G4VEnergyLossProcess, and G4VMultipleScattering.

Definition at line 124 of file G4VContinuousDiscreteProcess.cc.

References G4VProcess::ClearNumberOfInteractionLengthLeft(), and G4VProcess::pParticleChange.

00128 { 
00129 //  clear  NumberOfInteractionLengthLeft
00130     ClearNumberOfInteractionLengthLeft();
00131     return pParticleChange;
00132 }

G4double G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength ( const G4Track ,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety,
G4GPILSelection selection 
) [virtual]

Implements G4VProcess.

Reimplemented in G4VEnergyLossProcess, and G4VMultipleScattering.

Definition at line 134 of file G4VContinuousDiscreteProcess.cc.

References CandidateForSelection, G4DynamicParticle::DumpInfo(), G4cout, G4endl, GetContinuousStepLimit(), G4Track::GetDynamicParticle(), G4Track::GetMaterial(), G4Material::GetName(), G4VProcess::GetProcessName(), and G4VProcess::verboseLevel.

00141 {
00142   // GPILSelection is set to defaule value of CandidateForSelection
00143   valueGPILSelection = CandidateForSelection;
00144 
00145   // get Step limit proposed by the process
00146   G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep, currentSafety);
00147 
00148   // set return value for G4GPILSelection
00149   *selection = valueGPILSelection;
00150 
00151 #ifdef G4VERBOSE
00152    if (verboseLevel>1){
00153     G4cout << "G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength ";
00154     G4cout << "[ " << GetProcessName() << "]" <<G4endl;
00155     track.GetDynamicParticle()->DumpInfo();
00156     G4cout << " in Material  " <<  track.GetMaterial()->GetName() <<G4endl;
00157     G4cout << "IntractionLength= " << steplength/cm <<"[cm] " <<G4endl;
00158   }
00159 #endif
00160   return  steplength ;
00161 }

virtual G4VParticleChange* G4VContinuousDiscreteProcess::AtRestDoIt ( const G4Track ,
const G4Step  
) [inline, virtual]

Implements G4VProcess.

Definition at line 99 of file G4VContinuousDiscreteProcess.hh.

00102                               {return 0;};

virtual G4double G4VContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
) [inline, virtual]

Implements G4VProcess.

Definition at line 93 of file G4VContinuousDiscreteProcess.hh.

00096                               { return -1.0; };

virtual G4double G4VContinuousDiscreteProcess::GetContinuousStepLimit ( const G4Track aTrack,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety 
) [protected, pure virtual]

Implemented in G4VeLowEnergyLoss, G4hImpactIonisation, G4VEnergyLoss, G4VEnergyLossProcess, and G4VMultipleScattering.

Referenced by AlongStepGetPhysicalInteractionLength().

G4GPILSelection G4VContinuousDiscreteProcess::GetGPILSelection (  )  const [inline, protected]

Definition at line 128 of file G4VContinuousDiscreteProcess.hh.

00128 {return valueGPILSelection;};

virtual G4double G4VContinuousDiscreteProcess::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
) [protected, pure virtual]

Implemented in G4VeLowEnergyLoss, G4hImpactIonisation, G4hRDEnergyLoss, G4ePolarizedIonisation, G4VEnergyLoss, G4VEnergyLossProcess, and G4VMultipleScattering.

Referenced by PostStepGetPhysicalInteractionLength().

G4VParticleChange * G4VContinuousDiscreteProcess::PostStepDoIt ( const G4Track ,
const G4Step  
) [virtual]

Implements G4VProcess.

Reimplemented in G4VeLowEnergyLoss, G4hImpactIonisation, G4hRDEnergyLoss, G4VEnergyLoss, G4VEnergyLossProcess, and G4VMultipleScattering.

Definition at line 114 of file G4VContinuousDiscreteProcess.cc.

References G4VProcess::ClearNumberOfInteractionLengthLeft(), and G4VProcess::pParticleChange.

Referenced by G4hImpactIonisation::PostStepDoIt().

00118 { 
00119 //  clear  NumberOfInteractionLengthLeft
00120     ClearNumberOfInteractionLengthLeft();
00121     return pParticleChange;
00122 }

G4double G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
) [virtual]

Implements G4VProcess.

Reimplemented in G4ePolarizedIonisation, G4VEnergyLossProcess, and G4VMultipleScattering.

Definition at line 72 of file G4VContinuousDiscreteProcess.cc.

References G4VProcess::currentInteractionLength, DBL_MAX, G4DynamicParticle::DumpInfo(), G4cout, G4endl, G4Track::GetDynamicParticle(), G4Track::GetMaterial(), GetMeanFreePath(), G4Material::GetName(), G4VProcess::GetProcessName(), NotForced, G4VProcess::ResetNumberOfInteractionLengthLeft(), G4VProcess::SubtractNumberOfInteractionLengthLeft(), G4VProcess::theNumberOfInteractionLengthLeft, and G4VProcess::verboseLevel.

00077 {
00078   if ( (previousStepSize <=0.0) || (theNumberOfInteractionLengthLeft<=0.0)) {
00079     // beggining of tracking (or just after DoIt of this process)
00080     ResetNumberOfInteractionLengthLeft();
00081   } else if ( previousStepSize > 0.0) {
00082     // subtract NumberOfInteractionLengthLeft 
00083     SubtractNumberOfInteractionLengthLeft(previousStepSize);
00084   } else {
00085     // zero step
00086     //  DO NOTHING
00087   }
00088 
00089   // condition is set to "Not Forced"
00090   *condition = NotForced;
00091 
00092   // get mean free path
00093   currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
00094 
00095   G4double value;
00096   if (currentInteractionLength <DBL_MAX) {
00097     value = theNumberOfInteractionLengthLeft * currentInteractionLength;
00098   } else {
00099     value = DBL_MAX;
00100   }
00101 #ifdef G4VERBOSE
00102   if (verboseLevel>1){
00103     G4cout << "G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength ";
00104     G4cout << "[ " << GetProcessName() << "]" <<G4endl;
00105     track.GetDynamicParticle()->DumpInfo();
00106     G4cout << " in Material  " <<  track.GetMaterial()->GetName() <<G4endl;
00107     G4cout << "InteractionLength= " << value/cm <<"[cm] " <<G4endl;
00108   }
00109 #endif
00110   return value;
00111 }

void G4VContinuousDiscreteProcess::SetGPILSelection ( G4GPILSelection  selection  )  [inline, protected]

Definition at line 125 of file G4VContinuousDiscreteProcess.hh.

00126     { valueGPILSelection = selection;};


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:39 2013 for Geant4 by  doxygen 1.4.7