[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 | // $Id: G4VEmModel.hh,v 1.48 2007/10/29 08:38:58 vnivanch Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-01-patch-02 $ |
---|
| 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // GEANT4 Class header file |
---|
| 32 | // |
---|
| 33 | // |
---|
| 34 | // File name: G4VEmModel |
---|
| 35 | // |
---|
| 36 | // Author: Vladimir Ivanchenko |
---|
| 37 | // |
---|
| 38 | // Creation date: 03.01.2002 |
---|
| 39 | // |
---|
| 40 | // Modifications: |
---|
| 41 | // |
---|
| 42 | // 23-12-02 V.Ivanchenko change interface before move to cut per region |
---|
| 43 | // 24-01-03 Cut per region (V.Ivanchenko) |
---|
| 44 | // 13-02-03 Add name (V.Ivanchenko) |
---|
| 45 | // 25-02-03 Add sample theta and displacement (V.Ivanchenko) |
---|
| 46 | // 23-07-03 Replace G4Material by G4MaterialCutCouple in dE/dx and CrossSection |
---|
| 47 | // calculation (V.Ivanchenko) |
---|
| 48 | // 01-03-04 L.Urban signature changed in SampleCosineTheta |
---|
| 49 | // 23-04-04 L.urban signature of SampleCosineTheta changed back |
---|
| 50 | // 17-11-04 Add method CrossSectionPerAtom (V.Ivanchenko) |
---|
| 51 | // 14-03-05 Reduce number of pure virtual methods and make inline part |
---|
| 52 | // separate (V.Ivanchenko) |
---|
| 53 | // 24-03-05 Remove IsInCharge and add G4VParticleChange in the constructor (VI) |
---|
| 54 | // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko) |
---|
| 55 | // 15-04-05 optimize internal interface for msc (V.Ivanchenko) |
---|
| 56 | // 08-05-05 A -> N (V.Ivanchenko) |
---|
| 57 | // 25-07-05 Move constructor and destructor to the body (V.Ivanchenko) |
---|
| 58 | // 02-02-06 ComputeCrossSectionPerAtom: default value A=0. (mma) |
---|
| 59 | // 06-02-06 add method ComputeMeanFreePath() (mma) |
---|
| 60 | // 07-03-06 Optimize msc methods (V.Ivanchenko) |
---|
| 61 | // 29-06-06 Add member currentElement and Get/Set methods (V.Ivanchenko) |
---|
| 62 | // 29-10-07 Added SampleScattering (V.Ivanchenko) |
---|
| 63 | // |
---|
| 64 | // Class Description: |
---|
| 65 | // |
---|
| 66 | // Abstract interface to energy loss models |
---|
| 67 | |
---|
| 68 | // ------------------------------------------------------------------- |
---|
| 69 | // |
---|
| 70 | |
---|
| 71 | #ifndef G4VEmModel_h |
---|
| 72 | #define G4VEmModel_h 1 |
---|
| 73 | |
---|
| 74 | #include "globals.hh" |
---|
| 75 | #include "G4DynamicParticle.hh" |
---|
| 76 | #include "G4ParticleDefinition.hh" |
---|
| 77 | #include "G4MaterialCutsCouple.hh" |
---|
| 78 | #include "G4Material.hh" |
---|
| 79 | #include "G4Element.hh" |
---|
| 80 | #include "G4ElementVector.hh" |
---|
| 81 | #include "G4DataVector.hh" |
---|
| 82 | #include "G4VEmFluctuationModel.hh" |
---|
| 83 | #include "Randomize.hh" |
---|
| 84 | |
---|
| 85 | class G4PhysicsTable; |
---|
| 86 | class G4Region; |
---|
| 87 | class G4VParticleChange; |
---|
| 88 | class G4Track; |
---|
| 89 | |
---|
| 90 | class G4VEmModel |
---|
| 91 | { |
---|
| 92 | |
---|
| 93 | public: |
---|
| 94 | |
---|
| 95 | G4VEmModel(const G4String& nam); |
---|
| 96 | |
---|
| 97 | virtual ~G4VEmModel(); |
---|
| 98 | |
---|
| 99 | //------------------------------------------------------------------------ |
---|
| 100 | // Virtual methods to be implemented for the concrete model |
---|
| 101 | //------------------------------------------------------------------------ |
---|
| 102 | |
---|
| 103 | virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&) = 0; |
---|
| 104 | |
---|
| 105 | |
---|
| 106 | virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*, |
---|
| 107 | const G4MaterialCutsCouple*, |
---|
| 108 | const G4DynamicParticle*, |
---|
| 109 | G4double tmin = 0.0, |
---|
| 110 | G4double tmax = DBL_MAX) = 0; |
---|
| 111 | |
---|
| 112 | //------------------------------------------------------------------------ |
---|
| 113 | // Methods with standard implementation; may be overwritten if needed |
---|
| 114 | //------------------------------------------------------------------------ |
---|
| 115 | |
---|
| 116 | virtual G4double MinEnergyCut(const G4ParticleDefinition*, |
---|
| 117 | const G4MaterialCutsCouple*); |
---|
| 118 | |
---|
| 119 | virtual G4double ComputeDEDX(const G4MaterialCutsCouple*, |
---|
| 120 | const G4ParticleDefinition*, |
---|
| 121 | G4double kineticEnergy, |
---|
| 122 | G4double cutEnergy = DBL_MAX); |
---|
| 123 | |
---|
| 124 | virtual G4double CrossSection(const G4MaterialCutsCouple*, |
---|
| 125 | const G4ParticleDefinition*, |
---|
| 126 | G4double kineticEnergy, |
---|
| 127 | G4double cutEnergy = 0.0, |
---|
| 128 | G4double maxEnergy = DBL_MAX); |
---|
| 129 | |
---|
| 130 | virtual G4double ComputeDEDXPerVolume(const G4Material*, |
---|
| 131 | const G4ParticleDefinition*, |
---|
| 132 | G4double kineticEnergy, |
---|
| 133 | G4double cutEnergy = DBL_MAX); |
---|
| 134 | |
---|
| 135 | |
---|
| 136 | virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*, |
---|
| 137 | G4double kinEnergy, |
---|
| 138 | G4double Z, |
---|
| 139 | G4double A = 0., |
---|
| 140 | G4double cutEnergy = 0.0, |
---|
| 141 | G4double maxEnergy = DBL_MAX); |
---|
| 142 | |
---|
| 143 | virtual G4double ComputeMeanFreePath(const G4ParticleDefinition*, |
---|
| 144 | G4double kineticEnergy, |
---|
| 145 | const G4Material*, |
---|
| 146 | G4double cutEnergy = 0.0, |
---|
| 147 | G4double maxEnergy = DBL_MAX); |
---|
| 148 | |
---|
| 149 | virtual G4double CrossSectionPerVolume(const G4Material*, |
---|
| 150 | const G4ParticleDefinition*, |
---|
| 151 | G4double kineticEnergy, |
---|
| 152 | G4double cutEnergy = 0.0, |
---|
| 153 | G4double maxEnergy = DBL_MAX); |
---|
| 154 | |
---|
| 155 | protected: |
---|
| 156 | |
---|
| 157 | virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition*, |
---|
| 158 | G4double kineticEnergy); |
---|
| 159 | |
---|
| 160 | //------------------------------------------------------------------------ |
---|
| 161 | // Methods for msc simulation which needs to be overwritten |
---|
| 162 | //------------------------------------------------------------------------ |
---|
| 163 | |
---|
| 164 | public: |
---|
| 165 | |
---|
| 166 | virtual void SampleScattering(const G4DynamicParticle*, |
---|
| 167 | G4double safety); |
---|
| 168 | |
---|
| 169 | virtual G4double ComputeTruePathLengthLimit(const G4Track& track, |
---|
| 170 | G4PhysicsTable* theLambdaTable, |
---|
| 171 | G4double currentMinimalStep); |
---|
| 172 | |
---|
| 173 | virtual G4double ComputeGeomPathLength(G4double truePathLength); |
---|
| 174 | |
---|
| 175 | virtual G4double ComputeTrueStepLength(G4double geomPathLength); |
---|
| 176 | |
---|
| 177 | virtual void DefineForRegion(const G4Region*); |
---|
| 178 | |
---|
| 179 | //------------------------------------------------------------------------ |
---|
| 180 | // Generic methods common to all models |
---|
| 181 | //------------------------------------------------------------------------ |
---|
| 182 | |
---|
| 183 | void SetParticleChange(G4VParticleChange*, G4VEmFluctuationModel*); |
---|
| 184 | |
---|
| 185 | G4VEmFluctuationModel* GetModelOfFluctuations(); |
---|
| 186 | |
---|
| 187 | G4double HighEnergyLimit(); |
---|
| 188 | |
---|
| 189 | G4double LowEnergyLimit(); |
---|
| 190 | |
---|
| 191 | void SetHighEnergyLimit(G4double); |
---|
| 192 | |
---|
| 193 | void SetLowEnergyLimit(G4double); |
---|
| 194 | |
---|
| 195 | G4double MaxSecondaryKinEnergy(const G4DynamicParticle* dynParticle); |
---|
| 196 | |
---|
| 197 | const G4Element* SelectRandomAtom(const G4Material*, |
---|
| 198 | const G4ParticleDefinition*, |
---|
| 199 | G4double kineticEnergy, |
---|
| 200 | G4double cutEnergy = 0.0, |
---|
| 201 | G4double maxEnergy = DBL_MAX); |
---|
| 202 | |
---|
| 203 | const G4String& GetName() const; |
---|
| 204 | |
---|
| 205 | protected: |
---|
| 206 | |
---|
| 207 | const G4Element* GetCurrentElement() const; |
---|
| 208 | |
---|
| 209 | void SetCurrentElement(const G4Element*); |
---|
| 210 | |
---|
| 211 | private: |
---|
| 212 | |
---|
| 213 | // hide assignment operator |
---|
| 214 | G4VEmModel & operator=(const G4VEmModel &right); |
---|
| 215 | G4VEmModel(const G4VEmModel&); |
---|
| 216 | |
---|
| 217 | G4double lowLimit; |
---|
| 218 | G4double highLimit; |
---|
| 219 | G4double xsec[40]; |
---|
| 220 | |
---|
| 221 | G4VEmFluctuationModel* fluc; |
---|
| 222 | |
---|
| 223 | const G4String name; |
---|
| 224 | const G4Element* currentElement; |
---|
| 225 | |
---|
| 226 | protected: |
---|
| 227 | |
---|
| 228 | G4VParticleChange* pParticleChange; |
---|
| 229 | }; |
---|
| 230 | |
---|
| 231 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 232 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 233 | |
---|
| 234 | inline G4double G4VEmModel::HighEnergyLimit() |
---|
| 235 | { |
---|
| 236 | return highLimit; |
---|
| 237 | } |
---|
| 238 | |
---|
| 239 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 240 | |
---|
| 241 | inline G4double G4VEmModel::LowEnergyLimit() |
---|
| 242 | { |
---|
| 243 | return lowLimit; |
---|
| 244 | } |
---|
| 245 | |
---|
| 246 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 247 | |
---|
| 248 | inline void G4VEmModel::SetHighEnergyLimit(G4double val) |
---|
| 249 | { |
---|
| 250 | highLimit = val; |
---|
| 251 | } |
---|
| 252 | |
---|
| 253 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 254 | |
---|
| 255 | inline void G4VEmModel::SetLowEnergyLimit(G4double val) |
---|
| 256 | { |
---|
| 257 | lowLimit = val; |
---|
| 258 | } |
---|
| 259 | |
---|
| 260 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 261 | |
---|
| 262 | inline void G4VEmModel::SetParticleChange(G4VParticleChange* p, |
---|
| 263 | G4VEmFluctuationModel* f = 0) |
---|
| 264 | { |
---|
| 265 | if(p && pParticleChange != p) pParticleChange = p; |
---|
| 266 | fluc = f; |
---|
| 267 | } |
---|
| 268 | |
---|
| 269 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 270 | |
---|
| 271 | |
---|
| 272 | inline G4VEmFluctuationModel* G4VEmModel::GetModelOfFluctuations() |
---|
| 273 | { |
---|
| 274 | return fluc; |
---|
| 275 | } |
---|
| 276 | |
---|
| 277 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 278 | |
---|
| 279 | inline G4double G4VEmModel::MinEnergyCut(const G4ParticleDefinition*, |
---|
| 280 | const G4MaterialCutsCouple*) |
---|
| 281 | { |
---|
| 282 | return 0.0; |
---|
| 283 | } |
---|
| 284 | |
---|
| 285 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 286 | |
---|
| 287 | inline G4double G4VEmModel::ComputeDEDXPerVolume( |
---|
| 288 | const G4Material*, |
---|
| 289 | const G4ParticleDefinition*, |
---|
| 290 | G4double, |
---|
| 291 | G4double) |
---|
| 292 | { |
---|
| 293 | return 0.0; |
---|
| 294 | } |
---|
| 295 | |
---|
| 296 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 297 | |
---|
| 298 | inline G4double G4VEmModel::ComputeDEDX(const G4MaterialCutsCouple* c, |
---|
| 299 | const G4ParticleDefinition* p, |
---|
| 300 | G4double kinEnergy, |
---|
| 301 | G4double cutEnergy) |
---|
| 302 | { |
---|
| 303 | return ComputeDEDXPerVolume(c->GetMaterial(),p,kinEnergy,cutEnergy); |
---|
| 304 | } |
---|
| 305 | |
---|
| 306 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 307 | |
---|
| 308 | inline G4double G4VEmModel::CrossSection(const G4MaterialCutsCouple* c, |
---|
| 309 | const G4ParticleDefinition* p, |
---|
| 310 | G4double kinEnergy, |
---|
| 311 | G4double cutEnergy, |
---|
| 312 | G4double maxEnergy) |
---|
| 313 | { |
---|
| 314 | return |
---|
| 315 | CrossSectionPerVolume(c->GetMaterial(),p,kinEnergy,cutEnergy,maxEnergy); |
---|
| 316 | } |
---|
| 317 | |
---|
| 318 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 319 | |
---|
| 320 | inline G4double G4VEmModel::ComputeCrossSectionPerAtom( |
---|
| 321 | const G4ParticleDefinition*, |
---|
| 322 | G4double, G4double, G4double, |
---|
| 323 | G4double, G4double) |
---|
| 324 | { |
---|
| 325 | return 0.0; |
---|
| 326 | } |
---|
| 327 | |
---|
| 328 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 329 | |
---|
| 330 | inline const G4Element* G4VEmModel::SelectRandomAtom( |
---|
| 331 | const G4Material* material, |
---|
| 332 | const G4ParticleDefinition* pd, |
---|
| 333 | G4double kinEnergy, |
---|
| 334 | G4double tcut, |
---|
| 335 | G4double tmax) |
---|
| 336 | { |
---|
| 337 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
| 338 | G4int nelm = material->GetNumberOfElements(); |
---|
| 339 | currentElement = (*theElementVector)[nelm-1]; |
---|
| 340 | if (nelm > 1) { |
---|
| 341 | G4double x = G4UniformRand()* |
---|
| 342 | CrossSectionPerVolume(material,pd,kinEnergy,tcut,tmax); |
---|
| 343 | for(G4int i=0; i<nelm; i++) { |
---|
| 344 | if (x <= xsec[i]) { |
---|
| 345 | currentElement = (*theElementVector)[i]; |
---|
| 346 | break; |
---|
| 347 | } |
---|
| 348 | } |
---|
| 349 | } |
---|
| 350 | return currentElement; |
---|
| 351 | } |
---|
| 352 | |
---|
| 353 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 354 | |
---|
| 355 | inline const G4Element* G4VEmModel::GetCurrentElement() const |
---|
| 356 | { |
---|
| 357 | return currentElement; |
---|
| 358 | } |
---|
| 359 | |
---|
| 360 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 361 | |
---|
| 362 | inline void G4VEmModel::SetCurrentElement(const G4Element* elm) |
---|
| 363 | { |
---|
| 364 | currentElement = elm; |
---|
| 365 | } |
---|
| 366 | |
---|
| 367 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 368 | |
---|
| 369 | inline G4double G4VEmModel::MaxSecondaryKinEnergy( |
---|
| 370 | const G4DynamicParticle* dynParticle) |
---|
| 371 | { |
---|
| 372 | return MaxSecondaryEnergy(dynParticle->GetDefinition(), |
---|
| 373 | dynParticle->GetKineticEnergy()); |
---|
| 374 | } |
---|
| 375 | |
---|
| 376 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 377 | |
---|
| 378 | inline G4double G4VEmModel::MaxSecondaryEnergy(const G4ParticleDefinition*, |
---|
| 379 | G4double kineticEnergy) |
---|
| 380 | { |
---|
| 381 | return kineticEnergy; |
---|
| 382 | } |
---|
| 383 | |
---|
| 384 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 385 | |
---|
| 386 | inline const G4String& G4VEmModel::GetName() const |
---|
| 387 | { |
---|
| 388 | return name; |
---|
| 389 | } |
---|
| 390 | |
---|
| 391 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 392 | // Methods for msc simulation |
---|
| 393 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 394 | |
---|
| 395 | inline void G4VEmModel::SampleScattering(const G4DynamicParticle*, G4double) |
---|
| 396 | {} |
---|
| 397 | |
---|
| 398 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 399 | |
---|
| 400 | inline G4double G4VEmModel::ComputeTruePathLengthLimit( |
---|
| 401 | const G4Track&, |
---|
| 402 | G4PhysicsTable*, |
---|
| 403 | G4double) |
---|
| 404 | { |
---|
| 405 | return DBL_MAX; |
---|
| 406 | } |
---|
| 407 | |
---|
| 408 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 409 | |
---|
| 410 | inline G4double G4VEmModel::ComputeGeomPathLength(G4double truePathLength) |
---|
| 411 | { |
---|
| 412 | return truePathLength; |
---|
| 413 | } |
---|
| 414 | |
---|
| 415 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 416 | |
---|
| 417 | inline G4double G4VEmModel::ComputeTrueStepLength(G4double geomPathLength) |
---|
| 418 | { |
---|
| 419 | return geomPathLength; |
---|
| 420 | } |
---|
| 421 | |
---|
| 422 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 423 | |
---|
| 424 | inline void G4VEmModel::DefineForRegion(const G4Region*) |
---|
| 425 | {} |
---|
| 426 | |
---|
| 427 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 428 | |
---|
| 429 | #endif |
---|
| 430 | |
---|