00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039 #include "G4VContinuousDiscreteProcess.hh"
00040 #include "G4SystemOfUnits.hh"
00041
00042 #include "G4Step.hh"
00043 #include "G4Track.hh"
00044 #include "G4MaterialTable.hh"
00045 #include "G4VParticleChange.hh"
00046
00047 G4VContinuousDiscreteProcess::G4VContinuousDiscreteProcess()
00048 :G4VProcess("No Name Discrete Process"),
00049 valueGPILSelection(CandidateForSelection)
00050 {
00051 G4Exception("G4VContinuousDiscreteProcess::G4VContinuousDiscreteProcess()","ProcMan102",
00052 JustWarning,"Default constructor is called");
00053 }
00054
00055 G4VContinuousDiscreteProcess::G4VContinuousDiscreteProcess(const G4String& aName , G4ProcessType aType)
00056 : G4VProcess(aName, aType),
00057 valueGPILSelection(CandidateForSelection)
00058 {
00059 enableAtRestDoIt = false;
00060 }
00061
00062 G4VContinuousDiscreteProcess::~G4VContinuousDiscreteProcess()
00063 {
00064 }
00065
00066 G4VContinuousDiscreteProcess::G4VContinuousDiscreteProcess(G4VContinuousDiscreteProcess& right)
00067 : G4VProcess(right),
00068 valueGPILSelection(right.valueGPILSelection)
00069 {
00070 }
00071
00072 G4double G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(
00073 const G4Track& track,
00074 G4double previousStepSize,
00075 G4ForceCondition* condition
00076 )
00077 {
00078 if ( (previousStepSize <=0.0) || (theNumberOfInteractionLengthLeft<=0.0)) {
00079
00080 ResetNumberOfInteractionLengthLeft();
00081 } else if ( previousStepSize > 0.0) {
00082
00083 SubtractNumberOfInteractionLengthLeft(previousStepSize);
00084 } else {
00085
00086
00087 }
00088
00089
00090 *condition = NotForced;
00091
00092
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 }
00112
00113
00114 G4VParticleChange* G4VContinuousDiscreteProcess::PostStepDoIt(
00115 const G4Track& ,
00116 const G4Step&
00117 )
00118 {
00119
00120 ClearNumberOfInteractionLengthLeft();
00121 return pParticleChange;
00122 }
00123
00124 G4VParticleChange* G4VContinuousDiscreteProcess::AlongStepDoIt(
00125 const G4Track& ,
00126 const G4Step&
00127 )
00128 {
00129
00130 ClearNumberOfInteractionLengthLeft();
00131 return pParticleChange;
00132 }
00133
00134 G4double G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(
00135 const G4Track& track,
00136 G4double previousStepSize,
00137 G4double currentMinimumStep,
00138 G4double& currentSafety,
00139 G4GPILSelection* selection
00140 )
00141 {
00142
00143 valueGPILSelection = CandidateForSelection;
00144
00145
00146 G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep, currentSafety);
00147
00148
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 }