#include <G4IVRestDiscreteProcess.hh>
Inheritance diagram for G4IVRestDiscreteProcess:
Definition at line 59 of file G4IVRestDiscreteProcess.hh.
G4IVRestDiscreteProcess::G4IVRestDiscreteProcess | ( | const G4String & | , | |
G4ProcessType | aType = fNotDefined | |||
) |
Definition at line 49 of file G4IVRestDiscreteProcess.cc.
References G4VProcess::enableAlongStepDoIt, and G4IVRestDiscreteProcess().
Referenced by G4IVRestDiscreteProcess().
00050 : G4VProcess(aName, aType), 00051 theNlambdaTable(0),theInverseNlambdaTable(0), 00052 BIGSTEP(1.e10) 00053 { 00054 enableAlongStepDoIt = false; 00055 }
G4IVRestDiscreteProcess::G4IVRestDiscreteProcess | ( | G4IVRestDiscreteProcess & | ) |
Definition at line 61 of file G4IVRestDiscreteProcess.cc.
References G4IVRestDiscreteProcess().
00062 : G4VProcess(right), 00063 theNlambdaTable(0),theInverseNlambdaTable(0), 00064 BIGSTEP(1.e10) 00065 { 00066 }
G4IVRestDiscreteProcess::~G4IVRestDiscreteProcess | ( | ) | [virtual] |
virtual G4VParticleChange* G4IVRestDiscreteProcess::AlongStepDoIt | ( | const G4Track & | , | |
const G4Step & | ||||
) | [inline, virtual] |
virtual G4double G4IVRestDiscreteProcess::AlongStepGetPhysicalInteractionLength | ( | const G4Track & | , | |
G4double | , | |||
G4double | , | |||
G4double & | , | |||
G4GPILSelection * | ||||
) | [inline, virtual] |
G4VParticleChange * G4IVRestDiscreteProcess::AtRestDoIt | ( | const G4Track & | , | |
const G4Step & | ||||
) | [inline, virtual] |
Implements G4VProcess.
Definition at line 187 of file G4IVRestDiscreteProcess.hh.
References G4VProcess::ClearNumberOfInteractionLengthLeft(), and G4VProcess::pParticleChange.
00191 { 00192 // clear NumberOfInteractionLengthLeft 00193 ClearNumberOfInteractionLengthLeft(); 00194 00195 return pParticleChange; 00196 }
G4double G4IVRestDiscreteProcess::AtRestGetPhysicalInteractionLength | ( | const G4Track & | , | |
G4ForceCondition * | ||||
) | [inline, virtual] |
Implements G4VProcess.
Definition at line 159 of file G4IVRestDiscreteProcess.hh.
References G4VProcess::currentInteractionLength, G4DynamicParticle::DumpInfo(), G4cout, G4endl, G4Track::GetDynamicParticle(), G4Track::GetMaterial(), GetMeanLifeTime(), G4Material::GetName(), G4VProcess::GetProcessName(), NotForced, ns, G4VProcess::ResetNumberOfInteractionLengthLeft(), G4VProcess::theNumberOfInteractionLengthLeft, and G4VProcess::verboseLevel.
00163 { 00164 // beggining of tracking 00165 ResetNumberOfInteractionLengthLeft(); 00166 00167 // condition is set to "Not Forced" 00168 *condition = NotForced; 00169 00170 // get mean life time 00171 currentInteractionLength = GetMeanLifeTime(track, condition); 00172 00173 #ifdef G4VERBOSE 00174 if ((currentInteractionLength <0.0) || (verboseLevel>2)){ 00175 G4cout << "G4IVRestDiscreteProcess::AtRestGetPhysicalInteractionLength "; 00176 G4cout << "[ " << GetProcessName() << "]" <<G4endl; 00177 track.GetDynamicParticle()->DumpInfo(); 00178 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl; 00179 G4cout << "MeanLifeTime = " << currentInteractionLength/CLHEP::ns << "[ns]" <<G4endl; 00180 } 00181 #endif 00182 00183 return theNumberOfInteractionLengthLeft * currentInteractionLength; 00184 }
virtual G4double G4IVRestDiscreteProcess::GetMeanLifeTime | ( | const G4Track & | aTrack, | |
G4ForceCondition * | condition | |||
) | [protected, pure virtual] |
Referenced by AtRestGetPhysicalInteractionLength().
G4VParticleChange * G4IVRestDiscreteProcess::PostStepDoIt | ( | const G4Track & | , | |
const G4Step & | ||||
) | [inline, virtual] |
Implements G4VProcess.
Definition at line 148 of file G4IVRestDiscreteProcess.hh.
References G4VProcess::ClearNumberOfInteractionLengthLeft(), and G4VProcess::pParticleChange.
00152 { 00153 // reset NumberOfInteractionLengthLeft 00154 ClearNumberOfInteractionLengthLeft(); 00155 00156 return pParticleChange; 00157 }
G4double G4IVRestDiscreteProcess::PostStepGetPhysicalInteractionLength | ( | const G4Track & | track, | |
G4double | previousStepSize, | |||
G4ForceCondition * | condition | |||
) | [virtual] |
Implements G4VProcess.
Definition at line 68 of file G4IVRestDiscreteProcess.cc.
References DBL_MAX.
void G4IVRestDiscreteProcess::SubtractNumberOfInteractionLengthLeft | ( | G4double | previousStepSize | ) | [inline, protected, virtual] |
const G4double G4IVRestDiscreteProcess::BIGSTEP [protected] |
Definition at line 133 of file G4IVRestDiscreteProcess.hh.
Definition at line 131 of file G4IVRestDiscreteProcess.hh.
G4PhysicsTable* G4IVRestDiscreteProcess::theNlambdaTable [protected] |
Definition at line 130 of file G4IVRestDiscreteProcess.hh.