#include <G4ContinuousGainOfEnergy.hh>
Inheritance diagram for G4ContinuousGainOfEnergy:
Public Member Functions | |
G4ContinuousGainOfEnergy (const G4String &name="EnergyGain", G4ProcessType type=fElectromagnetic) | |
virtual | ~G4ContinuousGainOfEnergy () |
void | PreparePhysicsTable (const G4ParticleDefinition &) |
void | BuildPhysicsTable (const G4ParticleDefinition &) |
G4VParticleChange * | AlongStepDoIt (const G4Track &, const G4Step &) |
void | SetLossFluctuations (G4bool val) |
void | SetIsIntegral (G4bool val) |
void | SetDirectEnergyLossProcess (G4VEnergyLossProcess *aProcess) |
void | SetDirectParticle (G4ParticleDefinition *p) |
Protected Member Functions | |
virtual G4double | GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double ¤tSafety) |
Definition at line 69 of file G4ContinuousGainOfEnergy.hh.
G4ContinuousGainOfEnergy::G4ContinuousGainOfEnergy | ( | const G4String & | name = "EnergyGain" , |
|
G4ProcessType | type = fElectromagnetic | |||
) |
Definition at line 44 of file G4ContinuousGainOfEnergy.cc.
00045 : G4VContinuousProcess(name, type) 00046 { 00047 00048 00049 linLossLimit=0.05; 00050 lossFluctuationArePossible =true; 00051 lossFluctuationFlag=true; 00052 is_integral = false; 00053 00054 //Will be properly set in SetDirectParticle() 00055 IsIon=false; 00056 massRatio =1.; 00057 chargeSqRatio=1.; 00058 preStepChargeSqRatio=1.; 00059 00060 //Some initialization 00061 currentCoupleIndex=9999999; 00062 currentCutInRange=0.; 00063 currentMaterialIndex=9999999; 00064 currentTcut=0.; 00065 preStepKinEnergy=0.; 00066 preStepRange=0.; 00067 preStepScaledKinEnergy=0.; 00068 00069 currentCouple=0; 00070 }
G4ContinuousGainOfEnergy::~G4ContinuousGainOfEnergy | ( | ) | [virtual] |
G4VParticleChange * G4ContinuousGainOfEnergy::AlongStepDoIt | ( | const G4Track & | , | |
const G4Step & | ||||
) | [virtual] |
Reimplemented from G4VContinuousProcess.
Definition at line 114 of file G4ContinuousGainOfEnergy.cc.
References G4VProcess::aParticleChange, G4VEmModel::CorrectionsAlongStep(), G4VEmModel::GetChargeSquareRatio(), G4VEnergyLossProcess::GetDEDX(), G4Track::GetDynamicParticle(), G4VEnergyLossProcess::GetKineticEnergy(), G4DynamicParticle::GetKineticEnergy(), G4VEmModel::GetModelOfFluctuations(), G4Step::GetPostStepPoint(), G4VEnergyLossProcess::GetRange(), G4Step::GetStepLength(), G4StepPoint::GetWeight(), G4ParticleChange::Initialize(), G4VEmModel::MaxSecondaryKinEnergy(), CLHEP::detail::n, G4ParticleChange::ProposeEnergy(), G4VParticleChange::ProposeParentWeight(), G4DynamicParticle::SetDefinition(), G4VEnergyLossProcess::SetDynamicMassCharge(), G4DynamicParticle::SetKineticEnergy(), and G4VParticleChange::SetParentWeightByProcess().
00116 { 00117 00118 //Caution in this method the step length should be the true step length 00119 // A problem is that this is compute by the multiple scattering that does not know the energy at the end of the adjoint step. This energy is used during the 00120 //Forward sim. Nothing we can really do against that at this time. This is inherent to the MS method 00121 // 00122 00123 00124 00125 aParticleChange.Initialize(track); 00126 00127 // Get the actual (true) Step length 00128 //---------------------------------- 00129 G4double length = step.GetStepLength(); 00130 G4double degain = 0.0; 00131 00132 00133 00134 // Compute this for weight change after continuous energy loss 00135 //------------------------------------------------------------- 00136 G4double DEDX_before = theDirectEnergyLossProcess->GetDEDX(preStepKinEnergy, currentCouple); 00137 00138 00139 00140 // For the fluctuation we generate a new dynamic particle with energy =preEnergy+egain 00141 // and then compute the fluctuation given in the direct case. 00142 //----------------------------------------------------------------------- 00143 G4DynamicParticle* dynParticle = new G4DynamicParticle(); 00144 *dynParticle = *(track.GetDynamicParticle()); 00145 dynParticle->SetDefinition(theDirectPartDef); 00146 G4double Tkin = dynParticle->GetKineticEnergy(); 00147 00148 00149 size_t n=1; 00150 if (is_integral ) n=10; 00151 n=1; 00152 G4double dlength= length/n; 00153 for (size_t i=0;i<n;i++) { 00154 if (Tkin != preStepKinEnergy && IsIon) { 00155 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,Tkin); 00156 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio); 00157 00158 } 00159 00160 G4double r = theDirectEnergyLossProcess->GetRange(Tkin, currentCouple); 00161 if( dlength <= linLossLimit * r ) { 00162 degain = DEDX_before*dlength; 00163 } 00164 else { 00165 G4double x = r + dlength; 00166 //degain = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple) - theDirectEnergyLossProcess->GetKineticEnergy(r,currentCouple); 00167 G4double E = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple); 00168 if (IsIon){ 00169 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,E); 00170 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio); 00171 G4double x1= theDirectEnergyLossProcess->GetRange(E, currentCouple); 00172 while (std::abs(x-x1)>0.01*x) { 00173 E = theDirectEnergyLossProcess->GetKineticEnergy(x,currentCouple); 00174 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,E); 00175 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio); 00176 x1= theDirectEnergyLossProcess->GetRange(E, currentCouple); 00177 00178 } 00179 } 00180 00181 degain=E-Tkin; 00182 00183 00184 00185 } 00186 //G4cout<<degain<<G4endl; 00187 G4double tmax = currentModel->MaxSecondaryKinEnergy(dynParticle); 00188 tmax = std::min(tmax,currentTcut); 00189 00190 00191 dynParticle->SetKineticEnergy(Tkin+degain); 00192 00193 // Corrections, which cannot be tabulated for ions 00194 //---------------------------------------- 00195 G4double esecdep=0;//not used in most models 00196 currentModel->CorrectionsAlongStep(currentCouple, dynParticle, degain,esecdep, dlength); 00197 00198 // Sample fluctuations 00199 //------------------- 00200 00201 00202 G4double deltaE =0.; 00203 if (lossFluctuationFlag ) { 00204 deltaE = currentModel->GetModelOfFluctuations()-> 00205 SampleFluctuations(currentMaterial,dynParticle,tmax,dlength,degain)-degain; 00206 } 00207 00208 G4double egain=degain+deltaE; 00209 if (egain <=0) egain=degain; 00210 Tkin+=egain; 00211 dynParticle->SetKineticEnergy(Tkin); 00212 } 00213 00214 00215 00216 00217 00218 delete dynParticle; 00219 00220 if (IsIon){ 00221 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,Tkin); 00222 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatio); 00223 00224 } 00225 00226 G4double DEDX_after = theDirectEnergyLossProcess->GetDEDX(Tkin, currentCouple); 00227 00228 00229 G4double weight_correction=DEDX_after/DEDX_before; 00230 00231 00232 aParticleChange.ProposeEnergy(Tkin); 00233 00234 00235 //Caution!!! 00236 // It is important to select the weight of the post_step_point 00237 // as the current weight and not the weight of the track, as t 00238 // the weight of the track is changed after having applied all 00239 // the along_step_do_it. 00240 00241 // G4double new_weight=weight_correction*track.GetWeight(); //old 00242 G4double new_weight=weight_correction*step.GetPostStepPoint()->GetWeight(); 00243 aParticleChange.SetParentWeightByProcess(false); 00244 aParticleChange.ProposeParentWeight(new_weight); 00245 00246 00247 return &aParticleChange; 00248 00249 }
void G4ContinuousGainOfEnergy::BuildPhysicsTable | ( | const G4ParticleDefinition & | ) | [virtual] |
G4double G4ContinuousGainOfEnergy::GetContinuousStepLimit | ( | const G4Track & | track, | |
G4double | previousStepSize, | |||
G4double | currentMinimumStep, | |||
G4double & | currentSafety | |||
) | [protected, virtual] |
Implements G4VContinuousProcess.
Definition at line 262 of file G4ContinuousGainOfEnergy.cc.
References DBL_MAX, G4VEmModel::GetChargeSquareRatio(), G4Track::GetKineticEnergy(), G4Track::GetMaterialCutsCouple(), G4VEnergyLossProcess::GetRange(), G4VEmModel::HighEnergyLimit(), G4VEnergyLossProcess::SelectModelForMaterial(), and G4VEnergyLossProcess::SetDynamicMassCharge().
00264 { 00265 G4double x = DBL_MAX; 00266 x=.1*mm; 00267 00268 00269 DefineMaterial(track.GetMaterialCutsCouple()); 00270 00271 preStepKinEnergy = track.GetKineticEnergy(); 00272 preStepScaledKinEnergy = track.GetKineticEnergy()*massRatio; 00273 currentModel = theDirectEnergyLossProcess->SelectModelForMaterial(preStepScaledKinEnergy,currentCoupleIndex); 00274 G4double emax_model=currentModel->HighEnergyLimit(); 00275 if (IsIon) { 00276 chargeSqRatio = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,preStepKinEnergy); 00277 preStepChargeSqRatio = chargeSqRatio; 00278 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,preStepChargeSqRatio); 00279 } 00280 00281 00282 G4double maxE =1.1*preStepKinEnergy; 00283 /*if (preStepKinEnergy< 0.05*MeV) maxE =2.*preStepKinEnergy; 00284 else if (preStepKinEnergy< 0.1*MeV) maxE =1.5*preStepKinEnergy; 00285 else if (preStepKinEnergy< 0.5*MeV) maxE =1.25*preStepKinEnergy;*/ 00286 00287 if (preStepKinEnergy < currentTcut) maxE = std::min(currentTcut,maxE); 00288 00289 maxE=std::min(emax_model*1.001,maxE); 00290 00291 preStepRange = theDirectEnergyLossProcess->GetRange(preStepKinEnergy, currentCouple); 00292 00293 if (IsIon) { 00294 G4double chargeSqRatioAtEmax = currentModel->GetChargeSquareRatio(theDirectPartDef,currentMaterial,maxE); 00295 theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,chargeSqRatioAtEmax); 00296 } 00297 00298 G4double r1 = theDirectEnergyLossProcess->GetRange(maxE, currentCouple); 00299 00300 if (IsIon) theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,preStepChargeSqRatio); 00301 00302 00303 00304 x=r1-preStepRange; 00305 x=std::max(r1-preStepRange,0.001*mm); 00306 00307 return x; 00308 00309 00310 }
void G4ContinuousGainOfEnergy::PreparePhysicsTable | ( | const G4ParticleDefinition & | ) | [virtual] |
void G4ContinuousGainOfEnergy::SetDirectEnergyLossProcess | ( | G4VEnergyLossProcess * | aProcess | ) | [inline] |
void G4ContinuousGainOfEnergy::SetDirectParticle | ( | G4ParticleDefinition * | p | ) |
Definition at line 98 of file G4ContinuousGainOfEnergy.cc.
References G4ParticleDefinition::GetParticleType(), G4ParticleDefinition::GetPDGCharge(), and G4ParticleDefinition::GetPDGMass().
00099 {theDirectPartDef=p; 00100 if (theDirectPartDef->GetParticleType()== "nucleus") { 00101 IsIon=true; 00102 massRatio = proton_mass_c2/theDirectPartDef->GetPDGMass(); 00103 G4double q=theDirectPartDef->GetPDGCharge(); 00104 chargeSqRatio=q*q; 00105 00106 00107 } 00108 00109 }
void G4ContinuousGainOfEnergy::SetIsIntegral | ( | G4bool | val | ) | [inline] |
void G4ContinuousGainOfEnergy::SetLossFluctuations | ( | G4bool | val | ) |
Definition at line 252 of file G4ContinuousGainOfEnergy.cc.
00253 { 00254 if(val && !lossFluctuationArePossible) return; 00255 lossFluctuationFlag = val; 00256 }