File indexing completed on 2025-01-18 09:58:24
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054 #ifndef G4HIMPACTIONISATION
0055 #define G4HIMPACTIONISATION 1
0056
0057 #include <map>
0058 #include <CLHEP/Units/PhysicalConstants.h>
0059
0060 #include "globals.hh"
0061 #include "G4hRDEnergyLoss.hh"
0062 #include "G4DataVector.hh"
0063 #include "G4AtomicDeexcitation.hh"
0064 #include "G4PixeCrossSectionHandler.hh"
0065
0066 class G4VLowEnergyModel;
0067 class G4VParticleChange;
0068 class G4ParticleDefinition;
0069 class G4PhysicsTable;
0070 class G4MaterialCutsCouple;
0071 class G4Track;
0072 class G4Step;
0073
0074 class G4hImpactIonisation : public G4hRDEnergyLoss
0075 {
0076 public:
0077
0078 G4hImpactIonisation(const G4String& processName = "hImpactIoni");
0079
0080
0081
0082 ~G4hImpactIonisation();
0083
0084
0085 G4bool IsApplicable(const G4ParticleDefinition&);
0086
0087
0088 void BuildPhysicsTable(const G4ParticleDefinition& aParticleType) ;
0089
0090
0091 G4double GetMeanFreePath(const G4Track& track,
0092 G4double previousStepSize,
0093 enum G4ForceCondition* condition );
0094
0095
0096 void PrintInfoDefinition() const;
0097
0098
0099 void SetHighEnergyForProtonParametrisation(G4double energy) {protonHighEnergy = energy;} ;
0100
0101
0102
0103
0104 void SetLowEnergyForProtonParametrisation(G4double energy) {protonLowEnergy = energy;} ;
0105
0106
0107
0108
0109 void SetHighEnergyForAntiProtonParametrisation(G4double energy) {antiprotonHighEnergy = energy;} ;
0110
0111
0112
0113
0114 void SetLowEnergyForAntiProtonParametrisation(G4double energy) {antiprotonLowEnergy = energy;} ;
0115
0116
0117
0118
0119 G4double GetContinuousStepLimit(const G4Track& track,
0120 G4double previousStepSize,
0121 G4double currentMinimumStep,
0122 G4double& currentSafety);
0123
0124
0125 void SetElectronicStoppingPowerModel(const G4ParticleDefinition* aParticle,
0126 const G4String& dedxTable);
0127
0128
0129
0130 void SetNuclearStoppingPowerModel(const G4String& dedxTable)
0131 {theNuclearTable = dedxTable; SetNuclearStoppingOn();};
0132
0133
0134
0135
0136
0137
0138 void SetNuclearStoppingOn() {nStopping = true;};
0139
0140
0141 void SetNuclearStoppingOff() {nStopping = false;};
0142
0143
0144 void SetBarkasOn() {theBarkas = true;};
0145
0146
0147 void SetBarkasOff() {theBarkas = false;};
0148
0149
0150 void SetPixe(const G4bool ) {pixeIsActive = true;};
0151
0152
0153 G4VParticleChange* AlongStepDoIt(const G4Track& trackData ,
0154 const G4Step& stepData ) ;
0155
0156
0157 G4VParticleChange* PostStepDoIt(const G4Track& track,
0158 const G4Step& Step ) ;
0159
0160
0161 G4double ComputeDEDX(const G4ParticleDefinition* aParticle,
0162 const G4MaterialCutsCouple* couple,
0163 G4double kineticEnergy);
0164
0165
0166 void SetCutForSecondaryPhotons(G4double cut);
0167
0168
0169 void SetCutForAugerElectrons(G4double cut);
0170
0171
0172 void ActivateAugerElectronProduction(G4bool val);
0173
0174
0175
0176 void SetPixeCrossSectionK(const G4String& name) { modelK = name; }
0177 void SetPixeCrossSectionL(const G4String& name) { modelL = name; }
0178 void SetPixeCrossSectionM(const G4String& name) { modelM = name; }
0179 void SetPixeProjectileMinEnergy(G4double energy) { eMinPixe = energy; }
0180 void SetPixeProjectileMaxEnergy(G4double energy) { eMaxPixe = energy; }
0181
0182 protected:
0183
0184 private:
0185
0186 void InitializeMe();
0187 void InitializeParametrisation();
0188 void BuildLossTable(const G4ParticleDefinition& aParticleType);
0189
0190 void BuildLambdaTable(const G4ParticleDefinition& aParticleType);
0191 void SetProtonElectronicStoppingPowerModel(const G4String& dedxTable)
0192 {protonTable = dedxTable ;};
0193
0194
0195 void SetAntiProtonElectronicStoppingPowerModel(const G4String& dedxTable)
0196 {antiprotonTable = dedxTable;};
0197
0198 G4double MicroscopicCrossSection(const G4ParticleDefinition& aParticleType,
0199 G4double kineticEnergy,
0200 G4double atomicNumber,
0201 G4double deltaCutInEnergy) const;
0202
0203 G4double GetConstraints(const G4DynamicParticle* particle,
0204 const G4MaterialCutsCouple* couple);
0205
0206
0207 G4double ProtonParametrisedDEDX(const G4MaterialCutsCouple* couple,
0208 G4double kineticEnergy) const;
0209
0210 G4double AntiProtonParametrisedDEDX(const G4MaterialCutsCouple* couple,
0211 G4double kineticEnergy) const;
0212
0213 G4double DeltaRaysEnergy(const G4MaterialCutsCouple* couple,
0214 G4double kineticEnergy,
0215 G4double particleMass) const;
0216
0217
0218
0219 G4double BarkasTerm(const G4Material* material,
0220 G4double kineticEnergy) const;
0221
0222
0223 G4double BlochTerm(const G4Material* material,
0224 G4double kineticEnergy,
0225 G4double cSquare) const;
0226
0227
0228 G4double ElectronicLossFluctuation(const G4DynamicParticle* particle,
0229 const G4MaterialCutsCouple* material,
0230 G4double meanLoss,
0231 G4double step) const;
0232
0233
0234
0235 G4hImpactIonisation & operator=(const G4hImpactIonisation &right);
0236 G4hImpactIonisation(const G4hImpactIonisation&);
0237
0238 private:
0239
0240 G4VLowEnergyModel* betheBlochModel;
0241 G4VLowEnergyModel* protonModel;
0242 G4VLowEnergyModel* antiprotonModel;
0243 G4VLowEnergyModel* theIonEffChargeModel;
0244 G4VLowEnergyModel* theNuclearStoppingModel;
0245 G4VLowEnergyModel* theIonChuFluctuationModel;
0246 G4VLowEnergyModel* theIonYangFluctuationModel;
0247
0248
0249
0250
0251 G4String protonTable;
0252 G4String antiprotonTable;
0253 G4String theNuclearTable;
0254
0255
0256 G4double protonLowEnergy;
0257 G4double protonHighEnergy;
0258 G4double antiprotonLowEnergy;
0259 G4double antiprotonHighEnergy;
0260
0261
0262 G4bool nStopping;
0263 G4bool theBarkas;
0264
0265 G4DataVector cutForDelta;
0266 G4DataVector cutForGamma;
0267 G4double minGammaEnergy;
0268 G4double minElectronEnergy;
0269 G4PhysicsTable* theMeanFreePathTable;
0270
0271 const G4double paramStepLimit;
0272
0273 G4double fdEdx;
0274 G4double fRangeNow ;
0275 G4double charge;
0276 G4double chargeSquare;
0277 G4double initialMass;
0278 G4double fBarkas;
0279
0280 G4PixeCrossSectionHandler* pixeCrossSectionHandler;
0281 G4AtomicDeexcitation atomicDeexcitation;
0282 G4String modelK;
0283 G4String modelL;
0284 G4String modelM;
0285 G4double eMinPixe;
0286 G4double eMaxPixe;
0287
0288 G4bool pixeIsActive;
0289
0290 };
0291
0292
0293 inline G4double G4hImpactIonisation::GetContinuousStepLimit(const G4Track& track,
0294 G4double,
0295 G4double currentMinimumStep,
0296 G4double&)
0297 {
0298 G4double step = GetConstraints(track.GetDynamicParticle(),track.GetMaterialCutsCouple()) ;
0299
0300
0301
0302
0303
0304 if ((step > 0.) && (step < currentMinimumStep)) currentMinimumStep = step ;
0305
0306 return step ;
0307 }
0308
0309
0310 inline G4bool G4hImpactIonisation::IsApplicable(const G4ParticleDefinition& particle)
0311 {
0312
0313
0314
0315 return (particle.GetPDGCharge() != 0.0 && particle.GetPDGMass() > CLHEP::proton_mass_c2*0.1);
0316 }
0317
0318 #endif
0319
0320
0321
0322
0323
0324
0325
0326