#include <G4VRestContinuousDiscreteProcess.hh>
Inheritance diagram for G4VRestContinuousDiscreteProcess:
Definition at line 55 of file G4VRestContinuousDiscreteProcess.hh.
G4VRestContinuousDiscreteProcess::G4VRestContinuousDiscreteProcess | ( | const G4String & | , | |
G4ProcessType | aType = fNotDefined | |||
) |
Definition at line 55 of file G4VRestContinuousDiscreteProcess.cc.
References G4VRestContinuousDiscreteProcess().
Referenced by G4VRestContinuousDiscreteProcess().
00056 : G4VProcess(aName, aType), 00057 valueGPILSelection(CandidateForSelection) 00058 { 00059 }
G4VRestContinuousDiscreteProcess::G4VRestContinuousDiscreteProcess | ( | G4VRestContinuousDiscreteProcess & | ) |
Definition at line 65 of file G4VRestContinuousDiscreteProcess.cc.
References G4VRestContinuousDiscreteProcess().
00066 : G4VProcess(right), 00067 valueGPILSelection(right.valueGPILSelection) 00068 { 00069 }
G4VRestContinuousDiscreteProcess::~G4VRestContinuousDiscreteProcess | ( | ) | [virtual] |
G4VParticleChange * G4VRestContinuousDiscreteProcess::AlongStepDoIt | ( | const G4Track & | , | |
const G4Step & | ||||
) | [virtual] |
Implements G4VProcess.
Definition at line 139 of file G4VRestContinuousDiscreteProcess.cc.
References G4VProcess::pParticleChange.
00143 { 00144 return pParticleChange; 00145 }
G4double G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength | ( | const G4Track & | track, | |
G4double | previousStepSize, | |||
G4double | currentMinimumStep, | |||
G4double & | currentSafety, | |||
G4GPILSelection * | selection | |||
) | [virtual] |
Implements G4VProcess.
Definition at line 110 of file G4VRestContinuousDiscreteProcess.cc.
References CandidateForSelection, G4DynamicParticle::DumpInfo(), G4cout, G4endl, GetContinuousStepLimit(), G4Track::GetDynamicParticle(), G4Track::GetMaterial(), G4Material::GetName(), G4VProcess::GetProcessName(), and G4VProcess::verboseLevel.
00117 { 00118 // GPILSelection is set to defaule value of CandidateForSelection 00119 valueGPILSelection = CandidateForSelection; 00120 00121 // get Step limit proposed by the process 00122 G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep, currentSafety); 00123 00124 // set return value for G4GPILSelection 00125 *selection = valueGPILSelection; 00126 00127 #ifdef G4VERBOSE 00128 if (verboseLevel>1){ 00129 G4cout << "G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength "; 00130 G4cout << "[ " << GetProcessName() << "]" <<G4endl; 00131 track.GetDynamicParticle()->DumpInfo(); 00132 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl; 00133 G4cout << "IntractionLength= " << steplength/cm <<"[cm] " <<G4endl; 00134 } 00135 #endif 00136 return steplength ; 00137 }
G4VParticleChange * G4VRestContinuousDiscreteProcess::AtRestDoIt | ( | const G4Track & | , | |
const G4Step & | ||||
) | [virtual] |
Implements G4VProcess.
Definition at line 99 of file G4VRestContinuousDiscreteProcess.cc.
References G4VProcess::ClearNumberOfInteractionLengthLeft(), and G4VProcess::pParticleChange.
00103 { 00104 // clear NumberOfInteractionLengthLeft 00105 ClearNumberOfInteractionLengthLeft(); 00106 00107 return pParticleChange; 00108 }
G4double G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength | ( | const G4Track & | , | |
G4ForceCondition * | ||||
) | [virtual] |
Implements G4VProcess.
Definition at line 71 of file G4VRestContinuousDiscreteProcess.cc.
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.
00075 { 00076 // beggining of tracking 00077 ResetNumberOfInteractionLengthLeft(); 00078 00079 // condition is set to "Not Forced" 00080 *condition = NotForced; 00081 00082 // get mean life time 00083 currentInteractionLength = GetMeanLifeTime(track, condition); 00084 00085 #ifdef G4VERBOSE 00086 if ((currentInteractionLength <0.0) || (verboseLevel>2)){ 00087 G4cout << "G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength "; 00088 G4cout << "[ " << GetProcessName() << "]" <<G4endl; 00089 track.GetDynamicParticle()->DumpInfo(); 00090 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl; 00091 G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]" <<G4endl; 00092 } 00093 #endif 00094 00095 return theNumberOfInteractionLengthLeft * currentInteractionLength; 00096 }
virtual G4double G4VRestContinuousDiscreteProcess::GetContinuousStepLimit | ( | const G4Track & | aTrack, | |
G4double | previousStepSize, | |||
G4double | currentMinimumStep, | |||
G4double & | currentSafety | |||
) | [protected, pure virtual] |
Referenced by AlongStepGetPhysicalInteractionLength().
G4GPILSelection G4VRestContinuousDiscreteProcess::GetGPILSelection | ( | ) | const [inline, protected] |
virtual G4double G4VRestContinuousDiscreteProcess::GetMeanFreePath | ( | const G4Track & | aTrack, | |
G4double | previousStepSize, | |||
G4ForceCondition * | condition | |||
) | [protected, pure virtual] |
Referenced by PostStepGetPhysicalInteractionLength().
virtual G4double G4VRestContinuousDiscreteProcess::GetMeanLifeTime | ( | const G4Track & | aTrack, | |
G4ForceCondition * | condition | |||
) | [protected, pure virtual] |
Referenced by AtRestGetPhysicalInteractionLength().
G4VParticleChange * G4VRestContinuousDiscreteProcess::PostStepDoIt | ( | const G4Track & | , | |
const G4Step & | ||||
) | [virtual] |
Implements G4VProcess.
Definition at line 189 of file G4VRestContinuousDiscreteProcess.cc.
References G4VProcess::ClearNumberOfInteractionLengthLeft(), and G4VProcess::pParticleChange.
00193 { 00194 // clear NumberOfInteractionLengthLeft 00195 ClearNumberOfInteractionLengthLeft(); 00196 00197 return pParticleChange; 00198 }
G4double G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength | ( | const G4Track & | track, | |
G4double | previousStepSize, | |||
G4ForceCondition * | condition | |||
) | [virtual] |
Implements G4VProcess.
Definition at line 147 of file G4VRestContinuousDiscreteProcess.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.
00152 { 00153 if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0)) { 00154 // beggining of tracking (or just after DoIt of this process) 00155 ResetNumberOfInteractionLengthLeft(); 00156 } else if ( previousStepSize > 0.0) { 00157 // subtract NumberOfInteractionLengthLeft 00158 SubtractNumberOfInteractionLengthLeft(previousStepSize); 00159 } else { 00160 // zero step 00161 // DO NOTHING 00162 } 00163 00164 // condition is set to "Not Forced" 00165 *condition = NotForced; 00166 00167 // get mean free path 00168 currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition); 00169 00170 00171 G4double value; 00172 if (currentInteractionLength <DBL_MAX) { 00173 value = theNumberOfInteractionLengthLeft * currentInteractionLength; 00174 } else { 00175 value = DBL_MAX; 00176 } 00177 #ifdef G4VERBOSE 00178 if (verboseLevel>1){ 00179 G4cout << "G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength "; 00180 G4cout << "[ " << GetProcessName() << "]" <<G4endl; 00181 track.GetDynamicParticle()->DumpInfo(); 00182 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl; 00183 G4cout << "InteractionLength= " << value/cm <<"[cm] " <<G4endl; 00184 } 00185 #endif 00186 return value; 00187 }
void G4VRestContinuousDiscreteProcess::SetGPILSelection | ( | G4GPILSelection | selection | ) | [inline, protected] |