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 #include "G4ContinuousGainOfEnergy.hh"
00030
00031 #include "G4PhysicalConstants.hh"
00032 #include "G4SystemOfUnits.hh"
00033 #include "G4Step.hh"
00034 #include "G4ParticleDefinition.hh"
00035 #include "G4VEmModel.hh"
00036 #include "G4VEmFluctuationModel.hh"
00037 #include "G4VParticleChange.hh"
00038 #include "G4UnitsTable.hh"
00039 #include "G4AdjointCSManager.hh"
00040 #include "G4LossTableManager.hh"
00041
00043
00044 G4ContinuousGainOfEnergy::G4ContinuousGainOfEnergy(const G4String& name,
00045 G4ProcessType type): G4VContinuousProcess(name, type)
00046 {
00047
00048
00049 linLossLimit=0.05;
00050 lossFluctuationArePossible =true;
00051 lossFluctuationFlag=true;
00052 is_integral = false;
00053
00054
00055 IsIon=false;
00056 massRatio =1.;
00057 chargeSqRatio=1.;
00058 preStepChargeSqRatio=1.;
00059
00060
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 }
00071
00073
00074 G4ContinuousGainOfEnergy::~G4ContinuousGainOfEnergy()
00075 {
00076
00077 }
00079
00080
00081 void G4ContinuousGainOfEnergy::PreparePhysicsTable(
00082 const G4ParticleDefinition& )
00083 {
00084
00085 ;
00086 }
00087
00089
00090
00091 void G4ContinuousGainOfEnergy::BuildPhysicsTable(const G4ParticleDefinition&)
00092 {
00093 ;
00094 }
00095
00097
00098 void G4ContinuousGainOfEnergy::SetDirectParticle(G4ParticleDefinition* p)
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 }
00110
00112
00113
00114 G4VParticleChange* G4ContinuousGainOfEnergy::AlongStepDoIt(const G4Track& track,
00115 const G4Step& step)
00116 {
00117
00118
00119
00120
00121
00122
00123
00124
00125 aParticleChange.Initialize(track);
00126
00127
00128
00129 G4double length = step.GetStepLength();
00130 G4double degain = 0.0;
00131
00132
00133
00134
00135
00136 G4double DEDX_before = theDirectEnergyLossProcess->GetDEDX(preStepKinEnergy, currentCouple);
00137
00138
00139
00140
00141
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
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
00187 G4double tmax = currentModel->MaxSecondaryKinEnergy(dynParticle);
00188 tmax = std::min(tmax,currentTcut);
00189
00190
00191 dynParticle->SetKineticEnergy(Tkin+degain);
00192
00193
00194
00195 G4double esecdep=0;
00196 currentModel->CorrectionsAlongStep(currentCouple, dynParticle, degain,esecdep, dlength);
00197
00198
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
00236
00237
00238
00239
00240
00241
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 }
00251
00252 void G4ContinuousGainOfEnergy::SetLossFluctuations(G4bool val)
00253 {
00254 if(val && !lossFluctuationArePossible) return;
00255 lossFluctuationFlag = val;
00256 }
00258
00259
00260
00261
00262 G4double G4ContinuousGainOfEnergy::GetContinuousStepLimit(const G4Track& track,
00263 G4double , G4double , G4double& )
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
00284
00285
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 }
00311 #include "G4EmCorrections.hh"
00313
00314
00315 void G4ContinuousGainOfEnergy::SetDynamicMassCharge(const G4Track& ,G4double energy)
00316 {
00317
00318 G4double ChargeSqRatio= G4LossTableManager::Instance()->EmCorrections()->EffectiveChargeSquareRatio(theDirectPartDef,currentMaterial,energy);
00319 if (theDirectEnergyLossProcess) theDirectEnergyLossProcess->SetDynamicMassCharge(massRatio,ChargeSqRatio);
00320 }