[819] | 1 | // |
---|
| 2 | // ******************************************************************** |
---|
| 3 | // * License and Disclaimer * |
---|
| 4 | // * * |
---|
| 5 | // * The Geant4 software is copyright of the Copyright Holders of * |
---|
| 6 | // * the Geant4 Collaboration. It is provided under the terms and * |
---|
| 7 | // * conditions of the Geant4 Software License, included in the file * |
---|
| 8 | // * LICENSE and available at http://cern.ch/geant4/license . These * |
---|
| 9 | // * include a list of copyright holders. * |
---|
| 10 | // * * |
---|
| 11 | // * Neither the authors of this software system, nor their employing * |
---|
| 12 | // * institutes,nor the agencies providing financial support for this * |
---|
| 13 | // * work make any representation or warranty, express or implied, * |
---|
| 14 | // * regarding this software system or assume any liability for its * |
---|
| 15 | // * use. Please see the license in the file LICENSE and URL above * |
---|
| 16 | // * for the full disclaimer and the limitation of liability. * |
---|
| 17 | // * * |
---|
| 18 | // * This code implementation is the result of the scientific and * |
---|
| 19 | // * technical work of the GEANT4 collaboration. * |
---|
| 20 | // * By using, copying, modifying or distributing the software (or * |
---|
| 21 | // * any work based on the software) you agree to acknowledge its * |
---|
| 22 | // * use in resulting scientific publications, and indicate your * |
---|
| 23 | // * acceptance of all terms of the Geant4 Software license. * |
---|
| 24 | // ******************************************************************** |
---|
| 25 | // |
---|
| 26 | // |
---|
| 27 | // ------------------------------------------------------------ |
---|
| 28 | // GEANT 4 class header file |
---|
| 29 | // |
---|
| 30 | // History: based on object model of |
---|
| 31 | // 2nd December 1995, G.Cosmo |
---|
| 32 | // ---------- G4hLowEnergyIonisation physics process ----- |
---|
| 33 | // by Vladimir Ivanchenko, 14 July 1999 |
---|
| 34 | // was made on the base of G4hIonisation class |
---|
| 35 | // developed by Laszlo Urban |
---|
| 36 | // ************************************************************ |
---|
| 37 | |
---|
| 38 | // ************************************************************ |
---|
| 39 | // 28 July 1999 V.Ivanchenko cleen up |
---|
| 40 | // 17 August 1999 G.Mancinelli implemented ICRU parametrization (protons) |
---|
| 41 | // 20 August 1999 G.Mancinelli implemented ICRU parametrization (alpha) |
---|
| 42 | // 31 August 1999 V.Ivanchenko update and cleen up |
---|
| 43 | // 23 May 2000 MG Pia Clean up for QAO model |
---|
| 44 | // 25 July 2000 V.Ivanchenko New design iteration |
---|
| 45 | // 09 August 2000 V.Ivanchenko Add GetContinuousStepLimit |
---|
| 46 | // 17 August 2000 V.Ivanchenko Add IonFluctuationModel |
---|
| 47 | // 23 Oct 2000 V.Ivanchenko Renew comments |
---|
| 48 | // 30 Oct 2001 V.Ivanchenko Add minGammaEnergy and minElectronEnergy |
---|
| 49 | // 07 Dec 2001 V.Ivanchenko Add SetFluorescence method |
---|
| 50 | // 26 Feb 2002 V.Ivanchenko Add initialMass for GenericIons |
---|
| 51 | // 21 Jan 2003 V.Ivanchenko Cut per region |
---|
| 52 | // ------------------------------------------------------------ |
---|
| 53 | |
---|
| 54 | // Class Description: |
---|
| 55 | // Ionisation process of charged hadrons and ions, including low energy |
---|
| 56 | // extensions |
---|
| 57 | // The physics model is described in CERN-OPEN-99-121 and CERN-OPEN-99-300. |
---|
| 58 | // The user may select parametrisation tables for electronic |
---|
| 59 | // stopping powers and nuclear stopping powers |
---|
| 60 | // The list of available tables: |
---|
| 61 | // Electronic stopping powers: "ICRU_49p" (default), "ICRU_49He", |
---|
| 62 | // "Ziegler1977p", "Ziegler1985p", |
---|
| 63 | // "Ziegler1977He" |
---|
| 64 | // Nuclear stopping powers: "ICRU_49" (default), "Ziegler1977", |
---|
| 65 | // "Ziegler1985" |
---|
| 66 | // Further documentation available from http://www.ge.infn.it/geant4/lowE |
---|
| 67 | // and in the Physics Reference Manual |
---|
| 68 | |
---|
| 69 | // ------------------------------------------------------------ |
---|
| 70 | |
---|
| 71 | #ifndef G4hLowEnergyIonisation_h |
---|
| 72 | #define G4hLowEnergyIonisation_h 1 |
---|
| 73 | |
---|
| 74 | #include "globals.hh" |
---|
| 75 | #include "G4hLowEnergyLoss.hh" |
---|
| 76 | #include "G4VLowEnergyModel.hh" |
---|
| 77 | #include "G4Track.hh" |
---|
| 78 | #include "G4Step.hh" |
---|
| 79 | #include "G4Electron.hh" |
---|
| 80 | #include "G4PhysicsLogVector.hh" |
---|
| 81 | #include "G4PhysicsLinearVector.hh" |
---|
| 82 | #include "G4hNuclearStoppingModel.hh" |
---|
| 83 | #include "G4hBetheBlochModel.hh" |
---|
| 84 | #include "G4hParametrisedLossModel.hh" |
---|
| 85 | #include "G4QAOLowEnergyLoss.hh" |
---|
| 86 | #include "G4hIonEffChargeSquare.hh" |
---|
| 87 | #include "G4IonChuFluctuationModel.hh" |
---|
| 88 | #include "G4IonYangFluctuationModel.hh" |
---|
| 89 | #include "G4AtomicDeexcitation.hh" |
---|
| 90 | #include "G4MaterialCutsCouple.hh" |
---|
| 91 | #include <map> |
---|
| 92 | |
---|
| 93 | class G4VEMDataSet; |
---|
| 94 | class G4ShellVacancy; |
---|
| 95 | class G4VhShellCrossSection; |
---|
| 96 | |
---|
| 97 | class G4hLowEnergyIonisation : public G4hLowEnergyLoss |
---|
| 98 | { |
---|
| 99 | public: // With description |
---|
| 100 | |
---|
| 101 | G4hLowEnergyIonisation(const G4String& processName = "hLowEIoni"); |
---|
| 102 | // The ionisation process for hadrons/ions to be include in the |
---|
| 103 | // UserPhysicsList |
---|
| 104 | |
---|
| 105 | ~G4hLowEnergyIonisation(); |
---|
| 106 | // Destructor |
---|
| 107 | |
---|
| 108 | G4bool IsApplicable(const G4ParticleDefinition&); |
---|
| 109 | // True for all charged hadrons/ions |
---|
| 110 | |
---|
| 111 | void BuildPhysicsTable(const G4ParticleDefinition& aParticleType) ; |
---|
| 112 | // Build physics table during initialisation |
---|
| 113 | |
---|
| 114 | G4double GetMeanFreePath(const G4Track& track, |
---|
| 115 | G4double previousStepSize, |
---|
| 116 | enum G4ForceCondition* condition ); |
---|
| 117 | // Return MeanFreePath until delta-electron production |
---|
| 118 | |
---|
| 119 | void PrintInfoDefinition() const; |
---|
| 120 | // Print out of the class parameters |
---|
| 121 | |
---|
| 122 | void SetHighEnergyForProtonParametrisation(G4double energy) |
---|
| 123 | {protonHighEnergy = energy;} ; |
---|
| 124 | // Definition of the boundary proton energy. For higher energies |
---|
| 125 | // Bethe-Bloch formula is used, for lower energies a parametrisation |
---|
| 126 | // of the energy losses is performed. Default is 2 MeV. |
---|
| 127 | |
---|
| 128 | void SetLowEnergyForProtonParametrisation(G4double energy) |
---|
| 129 | {protonLowEnergy = energy;} ; |
---|
| 130 | // Set of the boundary proton energy. For lower energies |
---|
| 131 | // the Free Electron Gas model is used for the energy losses. |
---|
| 132 | // Default is 1 keV. |
---|
| 133 | |
---|
| 134 | void SetHighEnergyForAntiProtonParametrisation(G4double energy) |
---|
| 135 | {antiProtonHighEnergy = energy;} ; |
---|
| 136 | // Set of the boundary antiproton energy. For higher energies |
---|
| 137 | // Bethe-Bloch formula is used, for lower energies parametrisation |
---|
| 138 | // of the energy losses is performed. Default is 2 MeV. |
---|
| 139 | |
---|
| 140 | void SetLowEnergyForAntiProtonParametrisation(G4double energy) |
---|
| 141 | {antiProtonLowEnergy = energy;} ; |
---|
| 142 | // Set of the boundary antiproton energy. For lower energies |
---|
| 143 | // the Free Electron Gas model is used for the energy losses. |
---|
| 144 | // Default is 1 keV. |
---|
| 145 | |
---|
| 146 | G4double GetContinuousStepLimit(const G4Track& track, |
---|
| 147 | G4double previousStepSize, |
---|
| 148 | G4double currentMinimumStep, |
---|
| 149 | G4double& currentSafety); |
---|
| 150 | // Calculation of the step limit due to ionisation losses |
---|
| 151 | |
---|
| 152 | void SetElectronicStoppingPowerModel(const G4ParticleDefinition* aParticle, |
---|
| 153 | const G4String& dedxTable); |
---|
| 154 | // This method defines the electron ionisation parametrisation method |
---|
| 155 | // via the name of the table. Default is "ICRU_49p". |
---|
| 156 | |
---|
| 157 | void SetNuclearStoppingPowerModel(const G4String& dedxTable) |
---|
| 158 | {theNuclearTable = dedxTable; SetNuclearStoppingOn();}; |
---|
| 159 | // This method defines the nuclear ionisation parametrisation method |
---|
| 160 | // via the name of the table. Default is "ICRU_49". |
---|
| 161 | |
---|
| 162 | void SetNuclearStoppingOn() {nStopping = true;}; |
---|
| 163 | // This method switch on calculation of the nuclear stopping power. |
---|
| 164 | |
---|
| 165 | void SetNuclearStoppingOff() {nStopping = false;}; |
---|
| 166 | // This method switch off calculation of the nuclear stopping power. |
---|
| 167 | |
---|
| 168 | void SetBarkasOn() {theBarkas = true;}; |
---|
| 169 | // This method switch on calculation of the Barkas and Bloch effects. |
---|
| 170 | |
---|
| 171 | void SetBarkasOff() {theBarkas = false;}; |
---|
| 172 | // This method switch off calculation of the Barkas and Bloch effects. |
---|
| 173 | |
---|
| 174 | void SetFluorescence(const G4bool val) {theFluo = val;}; |
---|
| 175 | // This method switch on/off simulation of the fluorescence of the media. |
---|
| 176 | |
---|
| 177 | G4VParticleChange* AlongStepDoIt(const G4Track& trackData , |
---|
| 178 | const G4Step& stepData ) ; |
---|
| 179 | // Function to determine total energy deposition on the step |
---|
| 180 | |
---|
| 181 | G4VParticleChange* PostStepDoIt(const G4Track& track, |
---|
| 182 | const G4Step& Step ) ; |
---|
| 183 | // Simulation of delta rays production. |
---|
| 184 | |
---|
| 185 | G4double ComputeDEDX(const G4ParticleDefinition* aParticle, |
---|
| 186 | const G4MaterialCutsCouple* couple, |
---|
| 187 | G4double kineticEnergy); |
---|
| 188 | // This method returns electronic dE/dx for protons or antiproton. |
---|
| 189 | |
---|
| 190 | void SetCutForSecondaryPhotons(G4double cut); |
---|
| 191 | // Set threshold energy for fluorescence |
---|
| 192 | |
---|
| 193 | void SetCutForAugerElectrons(G4double cut); |
---|
| 194 | // Set threshold energy for Auger electron production |
---|
| 195 | |
---|
| 196 | void ActivateAugerElectronProduction(G4bool val); |
---|
| 197 | // Set Auger electron production flag on/off |
---|
| 198 | |
---|
| 199 | |
---|
| 200 | protected: |
---|
| 201 | |
---|
| 202 | private: |
---|
| 203 | |
---|
| 204 | void InitializeMe(); |
---|
| 205 | |
---|
| 206 | void InitializeParametrisation(); |
---|
| 207 | |
---|
| 208 | void BuildLossTable(const G4ParticleDefinition& aParticleType); |
---|
| 209 | |
---|
| 210 | void BuildDataForFluorescence(const G4ParticleDefinition& aParticleType); |
---|
| 211 | |
---|
| 212 | void BuildLambdaTable(const G4ParticleDefinition& aParticleType); |
---|
| 213 | |
---|
| 214 | void SetProtonElectronicStoppingPowerModel(const G4String& dedxTable) |
---|
| 215 | {theProtonTable = dedxTable ;}; |
---|
| 216 | // This method defines the ionisation parametrisation method via its name |
---|
| 217 | |
---|
| 218 | void SetAntiProtonElectronicStoppingPowerModel(const G4String& dedxTable) |
---|
| 219 | {theAntiProtonTable = dedxTable ;}; |
---|
| 220 | |
---|
| 221 | G4double ComputeMicroscopicCrossSection( |
---|
| 222 | const G4ParticleDefinition& aParticleType, |
---|
| 223 | G4double kineticEnergy, |
---|
| 224 | G4double atomicNumber, |
---|
| 225 | G4double deltaCutInEnergy) const; |
---|
| 226 | |
---|
| 227 | G4double GetConstraints(const G4DynamicParticle* particle, |
---|
| 228 | const G4MaterialCutsCouple* couple); |
---|
| 229 | // Function to determine StepLimit |
---|
| 230 | |
---|
| 231 | G4double ProtonParametrisedDEDX(const G4MaterialCutsCouple* couple, |
---|
| 232 | G4double kineticEnergy) const; |
---|
| 233 | |
---|
| 234 | G4double AntiProtonParametrisedDEDX(const G4MaterialCutsCouple* couple, |
---|
| 235 | G4double kineticEnergy) const; |
---|
| 236 | |
---|
| 237 | G4double DeltaRaysEnergy(const G4MaterialCutsCouple* couple, |
---|
| 238 | G4double kineticEnergy, |
---|
| 239 | G4double particleMass) const; |
---|
| 240 | // This method returns average energy loss due to delta-rays emission with |
---|
| 241 | // energy higher than the cut energy for given material. |
---|
| 242 | |
---|
| 243 | G4double BarkasTerm(const G4Material* material, |
---|
| 244 | G4double kineticEnergy) const; |
---|
| 245 | // Function to compute the Barkas term for protons |
---|
| 246 | |
---|
| 247 | G4double BlochTerm(const G4Material* material, |
---|
| 248 | G4double kineticEnergy, |
---|
| 249 | G4double cSquare) const; |
---|
| 250 | // Function to compute the Bloch term for protons |
---|
| 251 | |
---|
| 252 | G4double ElectronicLossFluctuation(const G4DynamicParticle* particle, |
---|
| 253 | const G4MaterialCutsCouple* material, |
---|
| 254 | G4double meanLoss, |
---|
| 255 | G4double step) const; |
---|
| 256 | // Function to sample electronic losses |
---|
| 257 | |
---|
| 258 | std::vector<G4DynamicParticle*>* DeexciteAtom(const G4MaterialCutsCouple* couple, |
---|
| 259 | G4double incidentEnergy, |
---|
| 260 | G4double hMass, |
---|
| 261 | G4double eLoss); |
---|
| 262 | |
---|
| 263 | G4int SelectRandomAtom(const G4MaterialCutsCouple* couple, |
---|
| 264 | G4double kineticEnergy) const; |
---|
| 265 | |
---|
| 266 | // hide assignment operator |
---|
| 267 | G4hLowEnergyIonisation & operator=(const G4hLowEnergyIonisation &right); |
---|
| 268 | G4hLowEnergyIonisation(const G4hLowEnergyIonisation&); |
---|
| 269 | |
---|
| 270 | private: |
---|
| 271 | // private data members ............................... |
---|
| 272 | G4VLowEnergyModel* theBetheBlochModel; |
---|
| 273 | G4VLowEnergyModel* theProtonModel; |
---|
| 274 | G4VLowEnergyModel* theAntiProtonModel; |
---|
| 275 | G4VLowEnergyModel* theIonEffChargeModel; |
---|
| 276 | G4VLowEnergyModel* theNuclearStoppingModel; |
---|
| 277 | G4VLowEnergyModel* theIonChuFluctuationModel; |
---|
| 278 | G4VLowEnergyModel* theIonYangFluctuationModel; |
---|
| 279 | std::map<G4int,G4double,std::less<G4int> > totalCrossSectionMap; |
---|
| 280 | |
---|
| 281 | // name of parametrisation table of electron stopping power |
---|
| 282 | G4String theProtonTable; |
---|
| 283 | G4String theAntiProtonTable; |
---|
| 284 | G4String theNuclearTable; |
---|
| 285 | |
---|
| 286 | // interval of parametrisation of electron stopping power |
---|
| 287 | G4double protonLowEnergy; |
---|
| 288 | G4double protonHighEnergy; |
---|
| 289 | G4double antiProtonLowEnergy; |
---|
| 290 | G4double antiProtonHighEnergy; |
---|
| 291 | |
---|
| 292 | // flag of parametrisation of nucleus stopping power |
---|
| 293 | G4bool nStopping; |
---|
| 294 | G4bool theBarkas; |
---|
| 295 | |
---|
| 296 | G4DataVector cutForDelta; |
---|
| 297 | G4DataVector cutForGamma; |
---|
| 298 | G4double minGammaEnergy; |
---|
| 299 | G4double minElectronEnergy; |
---|
| 300 | G4PhysicsTable* theMeanFreePathTable; |
---|
| 301 | |
---|
| 302 | const G4double paramStepLimit; // parameter limits the step at low energy |
---|
| 303 | |
---|
| 304 | G4double fdEdx; // computed in GetContraints |
---|
| 305 | G4double fRangeNow ; // |
---|
| 306 | G4double charge; // |
---|
| 307 | G4double chargeSquare; // |
---|
| 308 | G4double initialMass; // mass to calculate Lambda tables |
---|
| 309 | G4double fBarkas; |
---|
| 310 | |
---|
| 311 | G4AtomicDeexcitation deexcitationManager; |
---|
| 312 | G4ShellVacancy* shellVacancy; |
---|
| 313 | G4VhShellCrossSection* shellCS; |
---|
| 314 | std::vector<G4VEMDataSet*> zFluoDataVector; |
---|
| 315 | G4bool theFluo; |
---|
[1228] | 316 | G4bool expFlag; |
---|
[819] | 317 | }; |
---|
| 318 | |
---|
| 319 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 320 | |
---|
| 321 | inline G4double G4hLowEnergyIonisation::GetContinuousStepLimit( |
---|
| 322 | const G4Track& track, |
---|
| 323 | G4double, |
---|
| 324 | G4double currentMinimumStep, |
---|
| 325 | G4double&) |
---|
| 326 | { |
---|
| 327 | G4double Step = |
---|
| 328 | GetConstraints(track.GetDynamicParticle(),track.GetMaterialCutsCouple()) ; |
---|
| 329 | |
---|
| 330 | if((Step>0.0)&&(Step<currentMinimumStep)) |
---|
| 331 | currentMinimumStep = Step ; |
---|
| 332 | |
---|
| 333 | return Step ; |
---|
| 334 | } |
---|
| 335 | |
---|
| 336 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 337 | |
---|
| 338 | inline G4bool G4hLowEnergyIonisation::IsApplicable( |
---|
| 339 | const G4ParticleDefinition& particle) |
---|
| 340 | { |
---|
| 341 | return(particle.GetPDGCharge() != 0.0 |
---|
| 342 | && particle.GetPDGMass() > proton_mass_c2*0.1); |
---|
| 343 | } |
---|
| 344 | |
---|
| 345 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 346 | |
---|
| 347 | #endif |
---|
| 348 | |
---|
| 349 | |
---|
| 350 | |
---|
| 351 | |
---|
| 352 | |
---|
| 353 | |
---|
| 354 | |
---|
| 355 | |
---|